forked from keileg/fvbiot
-
Notifications
You must be signed in to change notification settings - Fork 0
/
mpsa.m
422 lines (329 loc) · 13.9 KB
/
mpsa.m
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
50
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
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
function out = mpsa(G,C,activeNodes,varargin)
% Discretize linear elasticity and Biot by a multipoint stress
% approximation.
%
% The implementation is modelled on mpfa.m, but everything becomes somewhat
% more complex for the vector unknown.
%
% NOTE: For poro-mechanical problems (not pure elasticity), non-homogeneous
% boundary conditions have not yet been implemented. This applies to the
% fields.
%
% Parameters:
% G - mrst grid structure, see http://www.sintef.no/projectweb/mrst/
% To get an overview of the data format, type 'help grid_structure'
% C - Constitutive relation. See shear_normal_stress.m for details.
% Active nodes: Which nodes should the stencil be computed for. If empty,
% all nodes in the grid is considered.
%
% NOTE: The function computes and stores in memory the inverse of a
% block-diagonal matrix
% (which expresses gradients on the sub-cells as a function of pressures
% in the surrounding cells). The size of each sub-matrix will be the
% number of sub-cells sharing the vertex * number of dimensions^2, so each
% system will be 16x16 for Cartesian 2D, 72x72 for Cartesian 3D. For large
% 3D grids, in particular with simplices, this may be heavy in terms of
% memory needs. For this reason, the activeNodes option may be an atractive
% option. For a work-around that splits the discretization into several
% calls to mpfa, and thus reduces memory consumption, see mpsa_subgrid.
%
% Optional parametrs, defined as keyword - value pairs:
% 'eta' - Location of continuity point on the half edges, defined
% according to Aavatsmark 2002. Between 0 and 1. Default value is 0, for
% simplices, the value should be 1/3 (will give symmetric method, see
% Klausen, Radu, Eigestad 2008).
% bc - boundary conditions, as defined by the MRST function addBC. If none
% are provided, homogeneous Neumann conditions will be assigned.
% invertBlocks - method for inverting block diagonal systems. Should be
% either 'matlab' or 'mex'. The former is pure matlab, which will be
% slow for large problems. Mex is substantially faster, but requires
% that matlab has access to a mex / c compiler.
% returnHalfEdgeValues - whether the return values for stresses etc should
% be for sub edges/faces (or summed to the whole face). Default is the
% full faces. Can be used for a partial update of the discretization
% scheme.
%
% Returns:
% out - a structure containing the following fields:
% Discretization of elasticity:
% A - discretization matrix
% div - divergence operator
% stress - discretized Hook's law
% boundStress - Discretization of boundary conditions
%
% Fields used for discretization of poro-mechanics
% divD - Discretized volumetric stress operator
% gradP - Discrete pressure gradient for fluid pressure in Biot's
% equations
% stabDelta - Stabilitazion term used for the flow equation.
%
% Example: See mpsa_ex and biot_ex
%
%{
Copyright 2015-2016, University of Bergen.
This file is part of FVBiot.
FVBiot is free software: you can redistribute it and/or modify
it under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
FVBiot 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 General Public License for more details.
You should have received a copy of the GNU General Public License
along with this file. If not, see <http://www.gnu.org/licenses/>.
%}
opt=struct('eta', 1/3,...
'bc',[],...
'invertBlocks','matlab',...
'verbose',false, ...
'returnHalfEdgeValues',0);
if nargin < 3 || numel(activeNodes) == 0
activeNodes = 1:G.nodes.num;
end
opt=merge_options(opt,varargin{:});
opt.invertBlocks = blockInverter(opt);
% Bookkeeping
Nd = G.griddim;
Nc = G.cells.num;
nFaceNodes = diff(G.faces.nodePos);
%% Boundary conditions
[isNeumann,isDirichlet] = classifyBoundaryFaces(G,opt.bc);
%% Subcell topology
[cnoAll, nnoAll, ~, fnoAll, subfnoAll, subhfnoAll] = createSubcellMapping(G);
localTopology = ismember(nnoAll, activeNodes);
cnoGlob = cnoAll(localTopology);
fnoGlob = fnoAll(localTopology);
subfnoGlob = subfnoAll(localTopology);
[~,~,nno] = unique(nnoAll(localTopology));
[tmpcno,~,cno] = unique(cnoAll(localTopology));
[~,~,fno] = unique(fnoAll(localTopology));
[~,~,subfno] = unique(subfnoAll(localTopology));
[~,~,subhfno] = unique(subhfnoAll(localTopology));
% Convert permeability tensor to matrix form for easy product with normal
% vector
C = C(tmpcno);
% Split the stiffness matrix into symmetric and asymmetric part.
[Csym,CAsym] = splitStiffnessMatrix(C);
clear C
%% Discretization
if opt.verbose
t = tic;
disp('Calculating nC products')
end
% Contribution of gradients to stress continuity (nC), both symmetric and
% asymmetric part
[nCsym,cellNodeBlocks,subcind] = multiplyNormalVectorStiffnessMatrix(G,Csym,cno,fnoGlob,nno,subhfno);
nCAsym = multiplyNormalVectorStiffnessMatrix(G,CAsym,cno,fnoGlob,nno,subhfno,1);
if opt.verbose
toc(t)
disp('done')
end
nsubhfno = numel(subhfno);
% Hook's law for normal vectors. This field is used only for Biot,
% specifically discretization of grad(p) (inclusion of fluid pressure in
% elasticity). For each face, it must be computed as seen from both
% neighboring cell, so do this before faces are uniquified (below)
[~,indf] = sort(reshape(repmat(subhfno,1,Nd),[],1));
hookNormal = sparse(1:(Nd*nsubhfno),indf,1) * (cat(1,nCsym{:}) + cat(1,nCAsym{:}));
nsubcells = size(cellNodeBlocks,1);
% Distance from cell centers to face continuity points
pCont = computeDistFaceCell(G,cnoGlob,fnoGlob,nno,subhfno,opt.eta);
% Unique subface numbers
[~,uniqueSubfno] = unique(subfno);
hfi = bsxfun(@minus,repmat(subfno(uniqueSubfno) * Nd,1,Nd),fliplr(1:Nd)-1);
fi = bsxfun(@minus,repmat(fno(uniqueSubfno) * Nd,1,Nd),fliplr(1:Nd)-1);
hf2f = sparse(fi,hfi,1);
hookSymCell = cell(Nd,1);
hookASymCell = hookSymCell;
for iter1 = 1 : Nd
hookSymCell{iter1} = nCsym{iter1}(uniqueSubfno,:);
hookASymCell{iter1} = nCAsym{iter1}(uniqueSubfno,:);
end
clear nCAsym
% nK and pCont was computed from both sides of a subface. Now combine them
sgn = 2*(cnoGlob == G.faces.neighbors(fnoGlob,1)) -1;
pairOverSubfaces = sparse(subfno,subhfno,sgn);
% Match stress equations from the two adjacent cells
for iter1 = 1 : Nd
nCsym{iter1} = pairOverSubfaces * nCsym{iter1};
end
pCont = pairOverSubfaces * pCont;
clear pairOverSubfaces
% Right hand side of the equations (corresponding to cell center
% contributions)
sgn = 2 * (cnoGlob == G.faces.neighbors(fnoGlob,1)) -1;
% Empty contribution from cell center displacements to face stresses (only
% strains count)
for iter1 = 1 : Nd
nKCC{iter1} = sparse(max(subfno),max(cno)*Nd);
end
pContCC = sparse(subfno,cno,sgn);
% Right hand side for fluid pressure gradient and stabilization term
rhsNormals = cell(Nd,1);
for iter1 = 1 : Nd
rhsNormals{iter1} = sparse(subfno,cno,G.faces.normals(fnoGlob,iter1) ...
.* sgn ./nFaceNodes(fnoGlob));
end
i = reshape(bsxfun(@minus,repmat(Nd*cno,1,Nd),fliplr(1:Nd)-1)',[],1);
j = (1:numel(subhfno)*Nd)';
v = reshape(repmat(sgn,1,Nd)',[],1);
divGradP = sparse(i,j,v);
clear i j v
% Done with matching expression from the two sides of the faces
nno = nno(uniqueSubfno);
subfno = subfno(uniqueSubfno);
nsubfno = max(subfno);
fnoGlob = fnoGlob(uniqueSubfno);
subfnoGlob = subfnoGlob(uniqueSubfno);
sgn = 2*(cnoGlob(uniqueSubfno) == G.faces.neighbors(fnoGlob,1)) -1;
% Now we can glue together Hook's law (symmetric and asymmetric parts) for
% the cell faces
hookCell = cell(Nd,1);
for iter1 = 1 : Nd
% Could have done this with cellfun
hookCell{iter1} = (hookSymCell{iter1} + hookASymCell{iter1});
end
% Final version of Hook's law. Multiply with gradients to get face stresses
[~,indf] = sort(repmat(subfno,Nd,1));
hook = sparse(1:Nd*numel(subfno),indf,1) * cat(1,hookCell{:});
clear hookSymCell hookASymCell hookCell indf
%% Exclude equations from boundary faces
j = find(~isNeumann(fnoGlob));
i = 1:numel(j);
excludeNeumann = sparse(i,j,1,numel(i),nsubfno);
j = find(~isDirichlet(fnoGlob));
i = 1 : numel(j);
excludeDirichlet = sparse(i,j,1,numel(i),nsubfno);
clear i j
for iter1 = 1 : Nd
nCsym{iter1} = excludeDirichlet * nCsym{iter1};
nKCC{iter1} = excludeDirichlet * nKCC{iter1};
end
pCont = excludeNeumann * pCont;
% Exclude boundary faces from right hand side
pContCC = excludeNeumann * pContCC;
%% Write the equations for the gradients on block diagonal form, and invert
% Map of rows to block-diag structure
% Do not consider excluded boundary equations
nnoStress = excludeDirichlet * nno;
nnoDispl = excludeNeumann * nno;
[n,map] = sort([repmat(nnoStress,Nd,1); repmat(nnoDispl,Nd,1)]);
rows2blkDiag = sparse(1:numel(map),map,1);
clear nnoStress nnoDispl map
% The columns in nC are number by covering all subcells in one cell at a
% time. To obtain a block structure, they must be centered around nodes
% instead
[~,sortNode] = sort(cellNodeBlocks(:,2));
subcind = reshape(subcind(sortNode,:)',[],1);
cols2blkDiag = sparse(subcind,1:numel(subcind),1);
nCsym = cat(1,nCsym{:}) * cols2blkDiag;
% Displacement contiuity conditions are originally formed for a scalar
% equation. First duplicate the conditions to get one for each component
pCont = kron(eye(Nd),pCont);
[~,map] = sort(repmat(reshape(ones(Nd,1) * (1:nsubcells),[],1),Nd,1));
pCont = pCont(:,map) * cols2blkDiag;
clear map
grad = rows2blkDiag * [nCsym ; pCont ];
nNeuCond = size(nCsym,1);
% Set up block matrix for gradient contribution
sz = accumarray(n,1);
% Clear some variables before we compute the inverse gradient
clear nCsym pCont
if opt.verbose
t = tic;
disp('Inverting gradient')
end
% Invert by chosen method
igrad = cols2blkDiag * opt.invertBlocks(grad,sz) * rows2blkDiag;
if opt.verbose
toc(t)
disp('done')
end
clear pCont grad sz n cols2blkDiag rows2blkDiag
%% Compute basis functions for internal degrees of freedom
% Rearrange rows of rhs, and scale the elements
% cc = rows2blkDiag * [nKCC *scalingNK; pContCC* scalingPcont];
nKCC = cat(1,nKCC{:});
pContCC = kron(eye(Nd),pContCC);
[~,map] = sort(repmat((1:numel(unique(cnoGlob)))',Nd,1));
pContCC = pContCC(:,map);
nDirCond = size(pContCC,1);
rif = bsxfun(@minus, Nd*repmat(subfnoGlob,1,Nd), fliplr(1:Nd)-1);
cif = bsxfun(@minus, Nd*repmat((1:numel(subfnoGlob))',1,Nd), fliplr(1:Nd)-1);
loc2globFace = sparse(rif,cif,1,Nd*max(subfnoAll),Nd*numel(subfnoGlob));
cnoGlob = unique(cnoGlob);
ric = bsxfun(@minus, Nd*repmat((1:numel(cnoGlob))',1,Nd), fliplr(1:Nd)-1);
cic = bsxfun(@minus, Nd*repmat(cnoGlob,1,Nd), fliplr(1:Nd)-1);
loc2globCell = sparse(ric, cic,1, Nd*numel(cnoGlob),Nd*max(cnoAll));
cc = [nKCC ; pContCC];
out.stress = loc2globFace * hook * igrad * (-cc) * loc2globCell;
clear pContCC nKCC
%% Boundary conditions
nNeuBnd = sum(isNeumann(fnoGlob)) * Nd;
nDirBnd = sum(isDirichlet(fnoGlob)) * Nd;
nbnd = nNeuBnd + nDirBnd;
fnoLoc = repmat(fnoGlob,Nd,1);
neuCol = find(excludeDirichlet * isNeumann(fnoGlob));
neuCol2 = neuCol;
neuCol = reshape(bsxfun(@plus,neuCol * ones(1,Nd),nNeuCond/Nd*((1:Nd)-1))',[],1);
neuVal = reshape(repmat(sgn(neuCol2),1,Nd)',[],1) .*G.faces.areas(fnoLoc(neuCol))...
./nFaceNodes(fnoLoc(neuCol));
ccNeu = sparse(neuCol,1:numel(neuCol),neuVal,nNeuCond,nbnd);
dirCol = find(excludeNeumann * isDirichlet(fnoGlob));
dirCol2 = dirCol;
dirCol = reshape(bsxfun(@plus,dirCol * ones(1,Nd),nDirCond/Nd*((1:Nd)-1))',[],1);
dirVal = reshape(repmat(sgn(dirCol2),1,Nd)',[],1);
ccDir = sparse(dirCol,nNeuBnd+(1:numel(dirCol)),dirVal,nDirCond,nNeuBnd+nDirBnd);
clear neuCol neuCol2 neuVal dirCol dirCol2 dirVal
tmp = [find(excludeDirichlet * isNeumann(fnoGlob)); find(excludeNeumann * isDirichlet(fnoGlob))];
cols = reshape(bsxfun(@minus,repmat(tmp * Nd,1,Nd),fliplr(1:Nd)-1)',[],1);
allFace2boundFace = sparse(1:nbnd,cols,1,nbnd,max(subfnoAll)*Nd);
map = @(i) reshape(bsxfun(@minus,repmat(i*Nd,1,Nd),fliplr(1:Nd)-1)',[],1);
face2halfFace = sparse(map(subfnoGlob),map(fnoGlob),1,Nd*max(subfnoAll),Nd*max(fnoAll));
out.boundStress = loc2globFace * hook * igrad * [ccNeu ; ccDir] * allFace2boundFace * face2halfFace;
clear hook face2halfFace allFace2boundFace tmp cols map
%% Coupling and stabilization terms for Biot
% Right hand side for the pressure gradient term
for iter1 = 1 : Nd
rhsNormals{iter1} = excludeDirichlet * rhsNormals{iter1};
end
rhsNormals = -[cat(1,rhsNormals{:}); sparse(nsubfno * Nd -nNeuBnd,numel(cnoGlob))];
loc2globCellsScalar = sparse(1:numel(cnoGlob),cnoGlob,1,numel(cnoGlob),Nc);
% Discretization of the term -grad(p) (impact of flow on elasticity)
out.gradP = loc2globCell' * divGradP * hookNormal * igrad * rhsNormals * ...
loc2globCellsScalar;
clear hookNormal
% Take divergence of mechanics and gradient to get divergence of
% displacement and stabilization term, respectively
if Nd == 2
trace = [1 4];
else
trace = [1 5 9];
end
[i,j] = ndgrid((1:nsubcells)',trace);
j = bsxfun(@plus,j,cumsum(Nd^2 * ones(nsubcells,1)) - Nd^2);
[~,cn] = gridCellNodes(G,cnoGlob);
cn = diff(cn);
cvol = G.cells.volumes(cnoGlob) ./ cn;
v = repmat(cvol(cellNodeBlocks(:,1)),1,Nd);
% Scale with subcell areas, and also add subcells together
divOp = sparse(cellNodeBlocks(:,1),(1:nsubcells)',1)* sparse(i,j,v);
% Impact of elasticity on flow equation (div(u))
out.divD = loc2globCellsScalar' * divOp * igrad * (-cc) * loc2globCell;
% Stabilization term, this is crucial when time step and compressibility
% goes to zero
out.stabDelta = -loc2globCellsScalar' * divOp * igrad * rhsNormals * loc2globCellsScalar;
clear igrad gradPsub i j cvol cn v rhsNormals
%%
% NOTE: Boundary conditions for the Biot terms have not been implemented.
% Should be straightforward extension of flow / elasticity.
%%
if opt.returnHalfEdgeValues
out.hf2f = hf2f;
else
out.stress = hf2f * out.stress;
out.boundStress = hf2f * out.boundStress;
out.div = vectorDivergence(G);
out.A = out.div * out.stress;
end