MATH661 Project 1
Large system factorization and model reduction

Posted: 10/02/23

Due: 10/16/21, 11:55PM (first draft, comments will be returned, revision due on 10/30)

This project presents a typical medium-size problem involving the concepts from approximation using linear combinations: transformation of coordinates, least squares, different types of basis sets.

1Introduction

1.1Modeling binary interaction in science and engineering

Many systems consist of point interactions. The iconic Eiffel Tower can be modeled through the trusses linking two points with coordinates 𝒙i,𝒙j3. The unit vector along direction from point i to point j is

𝒍ij=𝒙j-𝒙i||𝒙j-𝒙i||.

Suppose that structure deformation changes the coordinates to 𝒙i+𝒖i, 𝒙j+𝒖j. The most important structural response is a force 𝒇i exerted on point i that can be approximated as

𝒇i=𝒍ij𝒍ijT(𝒖j-𝒖i),

in appropriate physical units. An opposite force 𝒇j=-𝒇i acts on point j by principle of action and reaction. Note the projector along the truss direction 𝑷ij=𝒍ij𝒍ijT, allowing computation of nodal forces as

𝒇i=𝑷ij𝒖j-𝑷ij𝒖i
𝒇j=-𝑷ij𝒖j+𝑷ij𝒖i

Adding forces on a point from all trusses leads to the linear relation

𝒇=𝑲𝒖 (1)

with 𝒇,𝒖m the forces and displacements, 𝑲m×m the stiffness matrix, d=3 the number of dimensions, N the number of points, and m=dN, the number of degrees of freedom of the structure. This is a simple example of a finite element method. The vectors 𝒇,𝒖 contain the forces, displacements at all nodes, e.g.,

𝒖=[ 𝒖1 𝒖2 𝒖m ],

with each 𝒖id. Assembly of the stiffness matrix 𝑲 is carried out by adding contributions from each truss

[ 𝒇i 𝒇j ]=[ -𝑷ij 𝑷ij 𝑷ij -𝑷ij ][ 𝒖i 𝒖j ].

Such point interactions arise in molecular or social interactions, cellular motility, and economic exchange much in the same form (1). The techniques in this homework are thus equally applicable to physical chemistry, sociology, biology or business management.

Figure 1. (Left) Eiffel Tower model. (Center) Model nodes. (Right) Model edges.

1.2Problem data

Problem setup is carried through the following steps.

Add required Julia packages.

Load the data.

The following are typical instructions: open a file, load data, close the file.

ddir = homedir()*"/courses/MATH661/data/EiffelTower";
using MAT

Load the model point coordinates, display the first 3 node coordinates

dfile = matopen(ddir*"/EiffelPoints.mat");
X=read(dfile,"Expression1");
close(dfile);
NX,ndims=size(X);
X[1:3,:]

[ 12.563704490661621 -12.563704490661621 37.89190673828125 -12.563704490661621 -12.563704490661621 37.89190673828125 -14.107048034667969 -14.107048034667969 37.89190673828125 ] (2)

Load the model edges (links between the nodes) as floating point numbers (Matlab .mat file format)

dfile = matopen(ddir*"/EiffelLines.mat");
flL=read(dfile,"Expression1");
close(dfile);

Display model size

NL,nseg=size(flL); [NX ndims NL nseg]

[ 26464 3 31463 2 ] (3)

The model contains NX=26464 nodes (points) each with d=3 coordinates and NL edges (links) between node pairs. Transform the links to integers and display the first 4 edges.

L=zeros(Int64,NL,nseg);
for l=1:NL
  L[l,:] = floor.(Int64,flL[l,:])
end
L[1:4,:]

[ 519 520 520 521 521 519 540 541 ] (4)

The above states the first edge is from node i=519 to node j=520.

Define the system matrix 𝑲.

Figure 2. Structure of 𝑲 at various magnification factors

The data contains NX points in d, d=3, and NL line segments defined by 2 endpoints. The total number of degrees of freedom is m=dNX. Of the possible NX27×108 linkages between points in the structure, only NL=31463 are present.

m=3*NX; [NX NX^2 NL m m^2]

[ 26464 700343296 31463 79392 6303089664 ] (5)

The matrix 𝑲 can be formed by loops over the trusses as shown in the algorithm below.

Algorithm 𝑲 assembly

for k=1 to NL

i=Lk,1; j=Lk,2

𝒍=𝑿j-𝑿i; 𝒍=𝒍/||𝒍||; 𝑷=𝒍𝒍T

p=d(i-1); q=d(j-1)

for id=1 to d

for jd=1 to d

𝑲[p+id,p+jd]=𝑲[p+id,p+jd]-𝑷[id,jd]

𝑲[p+id,q+jd]=𝑲[p+id,q+jd]+𝑷[id,jd]

𝑲[q+id,p+jd]=𝑲[q+id,p+jd]+𝑷[id,jd]

𝑲[q+id,q+jd]=𝑲[q+id,q+jd]-𝑷[id,jd]

end

end

end

[ 𝒇i 𝒇j ]=[ -𝑷ij 𝑷ij 𝑷ij -𝑷ij ][ 𝒖i 𝒖j ].

Define a function to carry out the assembly process in d dimensions that returns a full matrix 𝑲.

function assembleK(d,L,X)
  NL = size(L)[1]; n=maximum(L); m=d*n; K=zeros(m,m)
  for k=1:NL
    i=L[k,1]; j=L[k,2]
    l=X[j,1:d]-X[i,1:d]; l=l/norm(l)
    p=d*(i-1); q=d*(j-1); P=l*l'
    for id=1:d
      for jd=1:d
        K[p+id,p+jd] = K[p+id,p+jd] - P[id,jd]
        K[p+id,q+jd] = K[p+id,q+jd] + P[id,jd]
        K[q+id,p+jd] = K[q+id,p+jd] + P[id,jd]
        K[q+id,q+jd] = K[q+id,q+jd] - P[id,jd]
      end
    end
  end
  return K
end;

The above has been precomputed for the Eiffel Tower model. It is wasteful, and often impossible due to memory constraints, to store the entire matrix of possible couplings between points 𝑲m×m (m26.3×109) since most of the matrix would consist of zero entries. Rather, only the nonzero elements corresponding to the linked nodes are stored through what is known as a sparse matrix. There are various techniques to store sparse matrices, one of which is the coordinate format defined by a vector of row indices i, a vector of column indices j and a vector of values v such that the (i,j) element of the matrix 𝑲 is

kij=K(i(k),j(k))=v(k).

The following instructions load the sparse representation of 𝑲.

dfile = matopen(ddir*"/KCOO.mat");
flrows=read(dfile,"rows");
flcols=read(dfile,"cols");
vals=read(dfile,"vals");
close(dfile);
nz=size(flrows)[1]

370621

rows=zeros(Int64,nz,1); cols=zeros(Int64,nz,1);
for l=1:nz
  rows[l] = floor(Int64,flrows[l])
  cols[l] = floor(Int64,flcols[l])
end
[minimum(rows) maximum(rows) minimum(cols) maximum(cols)]

[ 1 79392 1 79392 ] (6)

using SparseArrays
K=sparse(rows[:,1],cols[:,1],vals[:,1]);
size(K)

[ 79392 79392 ] (7)

K[5,6]

0.0

[rows[56] cols[56] vals[56]]

[ 22.0 14.0 1.9999999999999996 ] (8)

K[22,14]

1.9999999999999996

1.3Utility routines

Draw the deformed structure in figure nf, given node coordinates 𝑿NX×d, displacements 𝒖m, m=dNX, skipping ks nodes (to reduce drawing time).

Draw the deformed structure in figure nf, given node coordinates 𝑿NX×d, displacements 𝑼NX×d, m=dNX, skipping ks nodes.

1.4Simplified “A” model

The full Eiffel Tower model is quite complex. It is convenient to build familiarity with modeling binary interactions using the much simpler two-dimensional (d=2) 8-node model shown in Fig. 3

Figure 3. “A” model of Eiffel tower

1.5Enforcing boundary conditions

When structures are attached to their environment, boundary conditions have to be enforced in the linear system 𝒇=𝑲𝒖. A common technique is set the diagonal element of 𝑲 for the node i subjected to boundary conditions to a large value B such that it effectively decouples from the system, and set fi=Bgi where gi is the desired boundary value.

2Track 1 and 2 common problems

  1. Construct 𝑿A,𝑳A,𝑲A (node coordinates, edges, stiffness matrix) for the “A” model. For this small model a sparse representation of 𝑲A is not required.

  2. Start your investigation of the system by finding the n=10 largest singular values and eigenvalues of 𝑲,𝑲A. Plot the first p=3 singular modes and eigenmodes.

  3. Find the linear combination 𝒔 of the first n=10 singular modes that best approximates a force 𝒇=(zcosθ,zsinθ,0) on the structure by solving the problem

    min𝒔n||𝒇-𝑼𝒔||

    Randomly choose θ. This is the response of the structure to wind along direction θ. Form and solve the analogous problem for the “A” model first.

  4. Find the linear combination 𝒔 of the first n=10 eigenmodes modes that best approximates a force 𝒇=(zcosθ,zsinθ,0) on the structures. Form and solve the analogous problem for the “A” model first.

  5. Replace singular modes or eigenmodes by orthogonal 𝑸, with 𝑸𝑹=[ 𝒇 𝑲𝒇 𝑲2𝒇 𝑲9𝒇 ] (known as Krylov modes). Find the linear combination 𝒘 of the first 10 Krylov modes that best approximates a force 𝒇=(zcosθ,zsinθ,0) on the structure. Form and solve the analogous problem for the “A” model first.

3Track 2 task

Carry out a bibliographic search for review papers on reduced order modeling. Choose a paper you find appealing and draft a two-paragraph summary. Compare the approach you find in the literature to problems 3 to 5 above.