Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
ElmerCSC
GitHub Repository: ElmerCSC/elmerfem
Path: blob/devel/elmerice/Solvers/Flotation.F90
3204 views
1
!/*****************************************************************************/
2
! *
3
! * Elmer/Ice, a glaciological add-on to Elmer
4
! * http://elmerice.elmerfem.org
5
! *
6
! *
7
! * This program is free software; you can redistribute it and/or
8
! * modify it under the terms of the GNU General Public License
9
! * as published by the Free Software Foundation; either version 2
10
! * of the License, or (at your option) any later version.
11
! *
12
! * This program is distributed in the hope that it will be useful,
13
! * but WITHOUT ANY WARRANTY; without even the implied warranty of
14
! * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15
! * GNU General Public License for more details.
16
! *
17
! * You should have received a copy of the GNU General Public License
18
! * along with this program (in file fem/GPL-2); if not, write to the
19
! * Free Software Foundation, Inc., 51 Franklin Street, Fifth Floor,
20
! * Boston, MA 02110-1301, USA.
21
! *
22
! *****************************************************************************/
23
! ******************************************************************************
24
! *
25
! * Authors: F. GILLET-CHAULET
26
! * Email: [email protected]
27
! * Web: http://elmerice.elmerfem.org
28
! *
29
! * Original Date: Nov. 2017
30
! *
31
! *****************************************************************************
32
!!! Apply the flotation criterion to update Zb and Zs from the ice thickness H
33
!!! Execute this solver on the bottom boundary where the thickness equation is solved
34
!!! IF the mesh is vertically extruded export Zs on the top surface.
35
!
36
! OUTPUT Variables:
37
! Zb
38
! Zs
39
! DZbDt (optional, if variable found and transient simulation)
40
! DZsDt (optional, if variable found and transient simulation)
41
!
42
! INPUT Variable:
43
! H
44
! bedrock (optional)
45
!
46
! PARAMETERS:
47
! Constants:
48
! zsea
49
! rhow
50
! Material:
51
! SSA Mean Density
52
!
53
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
54
SUBROUTINE Flotation_init( Model,Solver,dt,TransientSimulation )
55
USE DefUtils
56
IMPLICIT NONE
57
!------------------------------------------------------------------------------
58
TYPE(Solver_t), TARGET :: Solver
59
TYPE(Model_t) :: Model
60
REAL(KIND=dp) :: dt
61
LOGICAL :: TransientSimulation
62
!--------------------------------------------------------------------------
63
CHARACTER(LEN=MAX_NAME_LEN) :: ZbName,ZsName,HName
64
CHARACTER(LEN=MAX_NAME_LEN) :: SolverName='Flotation'
65
TYPE(ValueList_t), POINTER :: SolverParams
66
LOGICAL :: GotIt
67
68
SolverParams => Solver % Values
69
70
CALL ListAddLogical( SolverParams,'No Matrix',.TRUE.)
71
72
ZbName = GetString(SolverParams, 'Bottom Surface Name', GotIt)
73
IF (.NOT.GotIt) THEN
74
CALL INFO(SolverName, 'Bottom Surface Name not found - using default Zb', level=3)
75
CALL ListAddString(SolverParams,'Bottom Surface Name','Zb')
76
END IF
77
78
ZsName = GetString(SolverParams, 'Top Surface Name', GotIt)
79
IF (.NOT.GotIt) THEN
80
CALL INFO(SolverName, 'Top Surface Name not found - using default Zs', level=3)
81
CALL ListAddString(SolverParams,'Top Surface Name','Zs')
82
END IF
83
84
HName = GetString(SolverParams, 'Thickness Variable Name', GotIt)
85
IF (.NOT.GotIt) THEN
86
CALL INFO(SolverName, 'Thickness Variable Name not found - using default H', level=3)
87
CALL ListAddString(SolverParams,'Thickness Variable Name','H')
88
END IF
89
90
END SUBROUTINE Flotation_init
91
!------------------------------------------------------------------------------
92
!------------------------------------------------------------------------------
93
SUBROUTINE Flotation( Model,Solver,dt,Transient )
94
!------------------------------------------------------------------------------
95
USE CoordinateSystems
96
USE MeshUtils
97
USE DefUtils
98
99
IMPLICIT NONE
100
!------------------------------------------------------------------------------
101
TYPE(Model_t) :: Model
102
TYPE(Solver_t), TARGET :: Solver
103
LOGICAL :: Transient
104
REAL(KIND=dp) :: dt
105
!------------------------------------------------------------------------------
106
! Local variables
107
!------------------------------------------------------------------------------
108
TYPE(Mesh_t),POINTER :: Mesh
109
TYPE(Solver_t),POINTER :: PSolver
110
TYPE(Variable_t),POINTER :: Var
111
TYPE(Variable_t),POINTER :: ZbVar,ZsVar
112
TYPE(Variable_t),POINTER :: HVar,BedVar
113
TYPE(Variable_t),POINTER :: GLMask,HafVar
114
TYPE(Variable_t),POINTER :: sftgif,sftgrf,sftflf
115
TYPE(Element_t), POINTER :: Element
116
TYPE(ValueList_t),POINTER :: BodyForce,Material, Params
117
TYPE(Nodes_t),SAVE :: ElementNodes
118
TYPE(GaussIntegrationPoints_t) :: IP
119
120
REAL(KIND=dp),DIMENSION(:),ALLOCATABLE,SAVE :: Density
121
REAL(KIND=dp),DIMENSION(:),ALLOCATABLE,SAVE :: GL
122
REAL(KIND=dp),DIMENSION(:),ALLOCATABLE,SAVE :: MinH,NodalH
123
REAL(KIND=dp),ALLOCATABLE,SAVE :: Basis(:)
124
REAL(KIND=dp) :: zsea,rhow,rhoi
125
REAL(KIND=dp) :: H,zb,zs,bedrock,Hf
126
REAL(KIND=dp), PARAMETER :: EPS=EPSILON(1.0)
127
REAL(KIND=dp) :: FillValue
128
REAL(KIND=dp) :: flarea,grarea,area,IParea,detJ
129
130
131
INTEGER,DIMENSION(:),POINTER,SAVE :: BotPointer,TopPointer,UpPointer
132
INTEGER, POINTER,SAVE :: NodeIndexes(:)
133
INTEGER :: GroundedNode
134
INTEGER :: ActiveDirection
135
INTEGER :: t,i,n,kk,ll
136
INTEGER :: topnode
137
INTEGER :: Active
138
INTEGER :: Eindex
139
INTEGER :: GlnIP
140
141
LOGICAL,SAVE :: Initialized = .FALSE.
142
LOGICAL,SAVE :: ExtrudedMesh=.False.
143
LOGICAL :: Found,GotIt
144
LOGICAL :: BoundarySolver
145
LOGICAL :: ComputeIceMasks,LimitedSolution,IceFree
146
LOGICAL :: stat
147
LOGICAL :: SEP
148
149
CHARACTER(LEN=MAX_NAME_LEN) :: SolverName='Flotation'
150
CHARACTER(LEN=MAX_NAME_LEN) :: ZbName,ZsName,HName
151
152
!------------------------------------------------------------------------------
153
154
Mesh => Model % Mesh
155
156
Params => Solver % Values
157
158
BoundarySolver = ( Solver % ActiveElements(1) > Model % Mesh % NumberOfBulkElements )
159
160
!!! get required variables Zb,Zs,H
161
ZbName = ListGetString(Params, 'Bottom Surface Name', UnFoundFatal=.TRUE.)
162
zbVar => VariableGet( Model % Mesh % Variables, ZbName,UnFoundFatal=.TRUE.)
163
164
ZsName = ListGetString(Params, 'Top Surface Name', UnFoundFatal=.TRUE.)
165
zsVar => VariableGet( Model % Mesh % Variables, ZsName,UnFoundFatal=.TRUE.)
166
167
HName = ListGetString(Params, 'Thickness Variable Name', UnFoundFatal=.TRUE.)
168
HVar => VariableGet( Model % Mesh % Variables, HName, UnFoundFatal=.TRUE.)
169
170
!!
171
!! get optional variables GLMAsk,bedrock
172
GLMAsk => VariableGet( Model % Mesh % Variables, 'GroundedMask')
173
IF (.NOT.ASSOCIATED(GLMAsk)) THEN
174
Message='GroundedMask not found'
175
CALL INFO(SolverName,Message,level=5)
176
ELSE
177
IF ((ParEnv % PEs>1).AND.(.NOT.ASSOCIATED(Solver%Matrix))) &
178
CALL FATAL(SolverName,'Solver%Matrix should be associated to update GLMask')
179
END IF
180
BedVar => VariableGet( Model % Mesh % Variables, 'bedrock')
181
IF (.NOT.ASSOCIATED(BedVar)) THEN
182
Message='bedrock not found'
183
CALL INFO(SolverName,Message,level=5)
184
END IF
185
! height above flotation
186
HafVar => VariableGet( Model % Mesh % Variables, 'Haf')
187
IF (.NOT.ASSOCIATED(HafVar)) THEN
188
Message='<Haf> not found; do not compute height above flotation'
189
CALL INFO(SolverName,Message,level=5)
190
END IF
191
192
!! compute ice area farctions
193
ComputeIceMasks = ListGetLogical(Params,"compute ice area fractions",Gotit)
194
LimitedSolution = .FALSE.
195
IF (ComputeIceMasks) THEN
196
FillValue = ListGetConstReal(Params,"Ice free mask values",Gotit)
197
IF (.NOT.Gotit) FillValue=0._dp
198
199
sftgif => VariableGet( Model % Mesh % Variables, 'sftgif',UnFoundFatal=.TRUE.)
200
sftgif % Values = 0._dp
201
IF (sftgif % TYPE /= Variable_on_elements) &
202
CALL FATAL(SolverName,"sftgif type should be on_elements")
203
sftgrf => VariableGet( Model % Mesh % Variables, 'sftgrf',UnFoundFatal=.TRUE.)
204
sftgrf % Values = 0._dp
205
IF (sftgrf % TYPE /= Variable_on_elements) &
206
CALL FATAL(SolverName,"sftgrf type should be on_elements")
207
sftflf => VariableGet( Model % Mesh % Variables, 'sftflf',UnFoundFatal=.TRUE.)
208
sftflf % Values = 0._dp
209
IF (sftflf % TYPE /= Variable_on_elements) &
210
CALL FATAL(SolverName,"sftflf type should be on_elements")
211
212
! change number of IPs for partially grounded elements
213
GLnIP=ListGetInteger( Params,'GL integration points number',SEP)
214
215
! check if we have a limited solution
216
! internal Elmer limiters...
217
LimitedSolution=ListCheckPresentAnyBodyForce(CurrentModel, TRIM(HName)//" Lower Limit")
218
IF (.NOT.LimitedSolution) &
219
LimitedSolution=ListCheckPresentAnyMaterial(CurrentModel, "Min "//TRIM(HName))
220
ENDIF
221
222
!!! Do some initialisation/allocation
223
IF ((.NOT.Initialized) .OR. Mesh % Changed) THEN
224
225
ActiveDirection = ListGetInteger(Params,'Active Coordinate',ExtrudedMesh)
226
IF (ExtrudedMesh) THEN
227
! Choose active direction coordinate and set corresponding unit vector
228
!---------------------------------------------------------------------
229
230
PSolver => Solver
231
CALL DetectExtrudedStructure( Mesh, PSolver, ExtVar = Var, BotNodePointer = BotPointer , &
232
TopNodePointer = TopPointer, UpNodePointer = UpPointer)
233
END IF
234
235
IF (Initialized) deallocate(Density,GL,MinH,NodalH,Basis)
236
237
N=Model % MaxElementNodes
238
allocate(Density(N),GL(N),MinH(N),NodalH(N),Basis(N))
239
240
Initialized = .TRUE.
241
END IF
242
!!
243
244
zsea = ListGetCReal( Model % Constants, 'Sea Level', UnFoundFatal=.TRUE. )
245
rhow = ListGetCReal( Model % Constants, 'water density', UnFoundFatal=.TRUE. )
246
247
IF (ASSOCIATED(GLMask)) GLMask%Values = -1.0
248
249
IF (BoundarySolver) THEN
250
Active = GetNOFBoundaryElements()
251
ELSE
252
Active = Solver % Mesh % NumberOfBulkElements
253
ENDIF
254
255
IF (ASSOCIATED(HafVar)) HafVar%Values = 0._dp
256
257
Do t=1,Active
258
259
IF (BoundarySolver) THEN
260
Element => GetBoundaryElement(t,Solver)
261
ELSE
262
Element => Solver % Mesh % Elements(t)
263
CurrentModel % CurrentElement => Element
264
ENDIF
265
266
Eindex = Element%ElementIndex
267
n = GetElementNOFNodes(Element)
268
NodeIndexes => Element % NodeIndexes
269
270
Material => GetMaterial(Element)
271
BodyForce => GetBodyForce(Element)
272
273
NodalH(1:n) = HVar%Values(HVar%Perm(NodeIndexes(1:n)))
274
275
IF (ComputeIceMasks) THEN
276
kk=sftgif % Perm(Eindex)
277
IF (kk>0) sftgif % Values ( kk ) = 1._dp
278
END IF
279
280
! If H is limited check if all element is at lower limit...
281
IceFree=.FALSE.
282
IF (LimitedSolution) THEN
283
MinH = 0._dp
284
! check for limited solution;
285
MinH = ListGetConstReal(BodyForce,TRIM(HName)//' Lower Limit', Gotit)
286
IF (.NOT.Gotit) &
287
MinH = ListGetConstReal(Material,'Min '//TRIM(HName), Gotit)
288
IF (.NOT.GotIt) CALL FATAL(SolverName,TRIM(HName)//" not found...but was supposed to be limited")
289
IF (ALL((NodalH(1:n)-MinH(1:n)) < EPS)) IceFree=.TRUE.
290
END IF
291
292
Density(1:n) = ListGetReal( Material, 'SSA Mean Density',n, NodeIndexes,UnFoundFatal=.TRUE.)
293
294
GroundedNode=0
295
GL=-1
296
Do i=1,n
297
298
H=NodalH(i)
299
rhoi=Density(i)
300
zb=zsea-H*rhoi/rhow
301
! if bedrock defined check flotation criterion
302
IF(ASSOCIATED(BedVar)) THEN
303
bedrock=BedVar%Values(BedVar%Perm(NodeIndexes(i)))
304
IF (zb < bedrock) THEN
305
zb=bedrock
306
GL(i)=1
307
GroundedNode=GroundedNode+1
308
END IF
309
IF (ASSOCIATED(HafVar)) THEN
310
Hf=MAX(zsea-bedrock,0._dp)*rhow/rhoi
311
HafVar%Values(HafVar%Perm(NodeIndexes(i)))=H-Hf
312
END IF
313
END IF
314
315
zs=zb+H
316
ZbVar%Values(ZbVar%Perm(NodeIndexes(i)))=zb
317
ZsVar%Values(ZsVar%Perm(NodeIndexes(i)))=zs
318
319
! Export Zs on top surface if required
320
IF (ExtrudedMesh) THEN
321
topnode=TopPointer(NodeIndexes(i))
322
ZsVar%Values(ZsVar%Perm(topnode))=zs
323
END IF
324
325
End do
326
IF (ASSOCIATED(GLMask)) THEN
327
IF ((GroundedNode > 0).AND.(GroundedNode < n)) THEN
328
WHERE(GL > 0._dp) GL=0._dp
329
END IF
330
GLMask%Values(GLMask%Perm(NodeIndexes(1:n)))= GL(1:n) * ABS(GLMask%Values(GLMask%Perm(NodeIndexes(1:n))))
331
END IF
332
333
334
IF (ComputeIceMasks) THEN
335
IF (GroundedNode == 0) THEN
336
!floating element
337
IF (sftgrf % Perm(Eindex) > 0 ) &
338
sftgrf % Values ( sftgrf % Perm(Eindex)) = 0._dp
339
IF (sftflf % Perm(Eindex) > 0 ) &
340
sftflf % Values ( sftflf % Perm(Eindex)) = 1._dp
341
ELSEIF (GroundedNode < n) THEN
342
!partly grounded
343
CALL GetElementNodes( ElementNodes, Element )
344
IF (SEP .AND. GLnIP > 0 ) THEN
345
IP = GaussPoints( Element ,np=GLnIP )
346
ELSE
347
IP = GaussPoints( Element )
348
ENDIF
349
area=0._dp
350
flarea=0._dp
351
grarea=0._dp
352
DO ll=1,IP % n
353
stat = ElementInfo( Element, ElementNodes, IP % U(ll), IP % V(ll), &
354
IP % W(ll), detJ, Basis )
355
rhoi=SUM(Density(1:n)*Basis(1:n))
356
H=SUM(NodalH(1:n)*Basis(1:n))
357
bedrock=SUM(BedVar%Values(BedVar%Perm(NodeIndexes(1:n)))*Basis(1:n))
358
Hf=(zsea-bedrock)*rhow/rhoi
359
IParea=detJ*IP % s(ll)
360
IF (H < Hf) THEN
361
flarea=flarea+IParea
362
ELSE
363
grarea=grarea+IParea
364
ENDIF
365
area=area+IParea
366
END DO
367
368
IF (sftgrf % Perm(Eindex) > 0 ) &
369
sftgrf % Values ( sftgrf % Perm(Eindex)) = grarea/area
370
IF (sftflf % Perm(Eindex) > 0 ) &
371
sftflf % Values ( sftflf % Perm(Eindex)) = flarea/area
372
ELSE
373
!grounded
374
IF (sftgrf % Perm(Eindex) > 0 ) &
375
sftgrf % Values ( sftgrf % Perm(Eindex)) = 1._dp
376
IF (sftflf % Perm(Eindex) > 0 ) &
377
sftflf % Values ( sftflf % Perm(Eindex)) = 0._dp
378
END IF
379
IF (IceFree) THEN
380
IF (sftgif % Perm(Eindex) > 0) &
381
sftgif % Values ( sftgif % Perm(Eindex)) = FillValue
382
IF (sftgrf % Perm(Eindex) > 0) &
383
sftgrf % Values ( sftgrf % Perm(Eindex)) = FillValue
384
IF (sftflf % Perm(Eindex) > 0) &
385
sftflf % Values ( sftflf % Perm(Eindex)) = FillValue
386
END IF
387
ENDIF
388
End Do
389
390
IF (ASSOCIATED(GLMask).AND.( ParEnv % PEs>1 )) CALL ParallelSumVector( Solver % Matrix, GLMask%Values ,1 )
391
392
!------------------------------------------------------------------------------
393
END SUBROUTINE Flotation
394
!------------------------------------------------------------------------------
395
396