Programming tasks to Scientific Computing I
Typedefs | Functions
ex1.cpp File Reference

Main file for 1st programming exercise. More...

#include <iostream>
#include <Eigen/Dense>
#include <Eigen/Sparse>
#include <Eigen/SparseCholesky>
#include <configuratorsShellFE.h>
#include <EigenDataContainer.h>
#include <quadratureShellFE.h>
#include <legacyVtkWriter.h>
#include <shellHandler.h>
#include <triangMesh.h>
#include <triangleShellFE.h>
#include "rhs.h"
#include "stiffnessMatrixIntegrator.h"

Go to the source code of this file.

Typedefs

typedef double RealType
 
typedef shellFE::ShellElementWithTangentSpaceAtVertex< DataTypeContainerShellFETriangleType
 
typedef shellFE::TriangMesh< DataTypeContainerShellFE, TriangleTypeMeshType
 
typedef shellFE::CenterOfMassQuadrature< RealType, typename DataTypeContainerShellFE::DomVecTypeQuadType
 
typedef shellFE::UnitTriangMeshConfiguratorP1< DataTypeContainerShellFE, MeshType, QuadTypeConfiguratorType
 
typedef ConfiguratorType::VectorType VectorType
 
typedef ConfiguratorType::SparseMatrixType SparseMatrixType
 

Functions

void printHelp ()
 
int main (int argc, char **argv)
 Main function containing the basic structure of the code. More...
 

Detailed Description

Main file for 1st programming exercise.

Definition in file ex1.cpp.

Typedef Documentation

Definition at line 26 of file ex1.cpp.

Definition at line 24 of file ex1.cpp.

Definition at line 25 of file ex1.cpp.

typedef double RealType

Definition at line 22 of file ex1.cpp.

Definition at line 28 of file ex1.cpp.

Definition at line 23 of file ex1.cpp.

Definition at line 27 of file ex1.cpp.

Function Documentation

int main ( int  argc,
char **  argv 
)

Main function containing the basic structure of the code.

Definition at line 39 of file ex1.cpp.

40 {
41  if (argc != 3) {
42  printHelp();
43  return 1;
44  }
45 
46  // Load a mesh
47  MeshType mesh(argv[1]);
48  // set neighbor relations
49  mesh.makeNeighbour();
50 
51  // The configurator handles information about basis functions,
52  // quadrature points etc.
53  ConfiguratorType conf(mesh);
54 
55  // The shellHandler handles additional information about the mesh
56  // including neighbor relations and the mapping from the unit triangle
57  // to the elements of the mesh
58  shellFE::ShellHandler<ConfiguratorType> shellHandler(conf);
59 
60  // Construct the right hand side
61  VectorType rhs(conf.getNumGlobalDofs());
62  createMaskedConstantRHS<ConfiguratorType>(conf, shellHandler, 1.0, rhs, true);
63 
64  // Construct the system matrix
65  SparseMatrixType systemMat(conf.getNumGlobalDofs(), conf.getNumGlobalDofs());
67  // Assemble the system matrix
68  integrand.assembleDirichlet<SparseMatrixType> (systemMat, shellHandler.getDirichletMask());
69 
70  // Construct a linear solver (Cholesky solver)
71  Eigen::SimplicialLLT<SparseMatrixType> solver;
72  // Decompose the matrix
73  solver.compute(systemMat);
74  // Solve the linear system
75  VectorType sol = solver.solve(rhs);
76 
77  // Save the mesh and solution as legacy VTK
78  LegacyVtkWriter<MeshType> vtkWriter ( mesh );
79  vtkWriter.addScalarData(sol, "solution", VERTEX_DATA);
80  vtkWriter.addScalarData(rhs, "rhs", VERTEX_DATA);
81  vtkWriter.save(argv[2]);
82 
83  return 0;
84 }
Additional information about TriangleMeshes.
Definition: shellHandler.h:28
void printHelp()
Definition: ex1.cpp:30
ConfiguratorType::SparseMatrixType SparseMatrixType
Definition: ex1.cpp:28
ConfiguratorType::VectorType VectorType
Definition: ex1.cpp:27
Configurator for Finite Elements.
void printHelp ( )

Definition at line 30 of file ex1.cpp.

31 {
32  cout << "Usage: ex1 meshFile resultFile\n";
33  cout << "\tmeshFile\tLegacy VTK file defining the Finite Element mesh "
34  << "to use\n";
35  cout << "\tresultFile\tLegacy VTK file to write the result to\n\n";
36 }