public static void TestDisplacements3x1()
        {
            // Expected solution
            double expectedStdUx5 = -8.17E-3;
            double expectedStdUx6 = -8.17E-3;

            double[] expectedEnrDisplNode5 = { 15.69E-3, 49.88E-3 };
            double[] expectedEnrDisplNode6 = { -15.69E-3, 49.88E-3 };

            // Analyze the model
            var dcb = new DoubleCantileverBeam();

            dcb.Create3x1Model();
            XModel model = dcb.Model;
            TrackingExteriorCrackLSM crack = dcb.Crack;

            (IVectorView globalU, IMatrixView globalK) = dcb.SolveModel();

            // Extract displacements of standard dofs
            DofTable freeDofs = model.Subdomains[DoubleCantileverBeam.subdomainID].FreeDofOrdering.FreeDofs;
            double   ux5      = globalU[freeDofs[model.Nodes[5], StructuralDof.TranslationX]];
            double   ux6      = globalU[freeDofs[model.Nodes[6], StructuralDof.TranslationX]];

            // Enriched dofs
            int numDofsPerNode = 2;
            IReadOnlyList <EnrichedDof> enrichedDofs = crack.CrackBodyEnrichment.Dofs;

            Assert.True(enrichedDofs.Count == numDofsPerNode);
            var enrDisplNode5 = new double[2];
            var enrDisplNode6 = new double[2];

            for (int i = 0; i < numDofsPerNode; ++i)
            {
                enrDisplNode5[i] = globalU[freeDofs[model.Nodes[5], enrichedDofs[i]]];
                enrDisplNode6[i] = globalU[freeDofs[model.Nodes[6], enrichedDofs[i]]];
            }

            // Check
            Func <double, double> round = x => 1E-3 * Math.Round(x * 1E3, 2);

            Assert.Equal(expectedStdUx5, round(ux5));
            Assert.Equal(expectedStdUx6, round(ux6));
            for (int i = 0; i < numDofsPerNode; ++i)
            {
                Assert.Equal(expectedEnrDisplNode5[i], round(enrDisplNode5[i]));
                Assert.Equal(expectedEnrDisplNode6[i], round(enrDisplNode6[i]));
            }
        }
        public static void TestStiffnesses3x1()
        {
            Matrix node6StiffnessExpected = 1E6 * Matrix.CreateFromArray(new double[, ]
            {
                { 1.154, 0.481, -0.481, -0.240 },
                { 0.481, 1.154, -0.240, -0.962 },
                { -0.481, -0.240, 0.962, 0.481 },
                { -0.240, -0.962, 0.481, 1.923 }
            });

            Matrix node7Elem1StiffnessExpected = 1E6 * Matrix.CreateFromArray(new double[, ]
            {
                { 1.154, 0.481, 1.568, 0.544, -0.444, 0.12, -0.847, 0.016, 0.378, -0.337 },
                { 0.481, 1.154, 0.575, 2.668, -0.114, 0.175, -0.165, -0.271, -0.055, -0.358 },
                { 1.568, 0.575, 12.432, 4.896, -0.824, -1.69, -3.114, -3.125, 0.134, 0.537 },
                { 0.544, 2.668, 4.896, 17.018, -1.359, -3.322, -2.444, -6.366, 0.465, 3.07 },
                { -0.444, -0.114, -0.824, -1.359, 2.639, -0.459, 1.648, -0.253, -1.869, 0.939 },
                { 0.12, 0.175, -1.69, -3.322, -0.459, 3.921, -0.214, 4.699, 0.909, -2.729 },
                { -0.847, -0.165, -3.114, -2.444, 1.648, -0.214, 4.063, -0.01, -1.386, 0.645 },
                { 0.016, -0.271, -3.125, -6.366, -0.253, 4.699, -0.01, 7.896, 0.685, -2.804 },
                { 0.378, -0.055, 0.134, 0.465, -1.869, 0.909, -1.386, 0.685, 3.081, -0.859 },
                { -0.337, -0.358, 0.537, 3.07, 0.939, -2.729, 0.645, -2.804, -0.859, 4.694 }
            });

            Matrix node7Elem2StiffnessExpected = 1E6 * Matrix.CreateFromArray(new double[, ]
            {
                { 1.154, -0.481, 1.715, -0.752, -0.586, 0.006, -0.886, -0.052, 0.600, -0.087 },
                { -0.481, 1.154, -0.856, 3.184, 0.092, -0.325, 0.106, -0.417, -0.156, 0.417 },
                { 1.715, -0.856, 13.871, -5.198, -2.050, 0.718, -3.537, 1.216, 1.561, -0.610 },
                { -0.752, 3.184, -5.198, 26.465, 0.679, -5.280, 1.135, -9.091, -0.753, 3.713 },
                { -0.586, 0.092, -2.050, 0.679, 1.098, -0.169, 1.822, -0.235, -0.843, 0.223 },
                { 0.006, -0.325, 0.718, -5.280, -0.169, 2.052, -0.246, 3.583, 0.234, -1.222 },
                { -0.886, 0.106, -3.537, 1.135, 1.822, -0.246, 3.075, -0.350, -1.317, 0.322 },
                { -0.052, -0.417, 1.216, -9.091, -0.235, 3.583, -0.350, 6.348, 0.330, -1.973 },
                { 0.600, -0.156, 1.561, -0.753, -0.843, 0.234, -1.317, 0.330, 0.869, -0.282 },
                { -0.087, 0.417, -0.610, 3.713, 0.223, -1.222, 0.322, -1.973, -0.282, 1.180 }
            });

            Matrix node7GlobalStiffnessExpected = 1E6 * Matrix.CreateFromArray(new double[, ]
            {
                { 2.308, 0.000, 3.283, -0.208, -1.030, 0.126, -1.733, -0.036, 0.979, -0.424 },
                { 0.000, 2.308, -0.282, 5.852, -0.022, -0.150, -0.060, -0.687, -0.211, 0.059 },
                { 3.283, -0.282, 26.303, -0.302, -2.874, -0.972, -6.651, -1.910, 1.695, -0.073 },
                { -0.208, 5.852, -0.302, 43.483, -0.680, -8.601, -1.310, -15.456, -0.289, 6.783 },
                { -1.030, -0.022, -2.874, -0.680, 3.736, -0.628, 3.470, -0.488, -2.712, 1.162 },
                { 0.126, -0.150, -0.972, -8.601, -0.628, 5.973, -0.460, 8.282, 1.142, -3.951 },
                { -1.733, -0.060, -6.651, -1.310, 3.470, -0.460, 7.139, -0.360, -2.703, 0.966 },
                { -0.036, -0.687, -1.910, -15.456, -0.488, 8.282, -0.360, 14.244, 1.015, -4.777 },
                { 0.979, -0.211, 1.695, -0.289, -2.712, 1.142, -2.703, 1.015, 3.950, -1.141 },
                { -0.424, 0.059, -0.073, 6.783, 1.162, -3.951, 0.966, -4.777, -1.141, 5.874 }
            });

            // Create and analyze model, in order to get the global stiffness
            var dcb = new DoubleCantileverBeam();

            dcb.Create3x1Model();
            XModel model = dcb.Model;
            TrackingExteriorCrackLSM crack = dcb.Crack;

            (IVectorView globalU, IMatrixView globalK) = dcb.SolveModel();

            // Print matrix
//            var writer = new FullMatrixWriter();
            //writer.NumericFormat = new FixedPointFormat() { NumDecimalDigits = 2 };
//            writer.ArrayFormat = new Array2DFormat("", "", "", "\n", ",");
//            writer.WriteToFile(globalK/*.DoToAllEntries(x => Math.Round(x * 1E-6, 3))*/, @"C:\Users\Serafeim\Desktop\xfem.txt");

            // Calculate relevant stiffness submatrix global
            XNode    node7           = model.Nodes[7];
            var      node7GlobalDofs = new int[10];
            DofTable freeDofs        = model.Subdomains[DoubleCantileverBeam.subdomainID].FreeDofOrdering.FreeDofs;

            node7GlobalDofs[0] = freeDofs[node7, StructuralDof.TranslationX];
            node7GlobalDofs[1] = freeDofs[node7, StructuralDof.TranslationY];
            for (int i = 0; i < 8; ++i)
            {
                node7GlobalDofs[2 + i] = freeDofs[node7, crack.CrackTipEnrichments.Dofs[i]];
            }
            IMatrix node7GlobalStiffness = globalK.GetSubmatrix(node7GlobalDofs, node7GlobalDofs);

            // Element 1 dofs (std first):
            // (N1,ux,0) (N1,uy,1) (N4,ux,2) (N4,uy,3) (N7,ux,4) (N7,uy,5) (N2,ux,6) (N2,uy,7)
            // (N1,tip0x, 8) (N1,tip0y,9) (N1,tip1x,10) (N1,tip1y,11) (N1,tip2x,12) (N1,tip2y,13) (N1,tip3x,14) (N1,tip3y,15)
            // (N4,tip0x,16) (N4,tip0y,17) (N4,tip1x,18) (N4,tip1y,19) (N4,tip2x,20) (N4,tip2y,21) (N4,tip3x,22) (N4,tip3y,23)
            // (N7,tip0x,24) (N7,tip0y,25) (N7,tip1x,26) (N7,tip1y,27) (N7,tip2x,28) (N7,tip2y,29) (N7,tip3x,30) (N7,tip3y,31)
            // (N2,tip0x,32) (N2,tip0y,33) (N2,tip1x,34) (N2,tip1y,35) (N2,tip2x,36) (N2,tip2y,37) (N2,tip3x,38) (N2,tip3y,39)

            // Element 2 dofs (std first):
            // (N4,ux,0) (N4,uy,1) (N5,ux,2) (N5,uy,3) (N6,ux,4) (N6,uy,5) (N7,ux,6) (N7,uy,7)
            // (N4,tip0x,8) (N4,tip0y,9) (N4,tip1x,10) (N4,tip1y,11) (N4,tip2x,12) (N4,tip2y,13) (N4,tip3x,14) (N4,tip3y,15)
            // (N5,bodyX,16) (N5,bodyY,17)
            // (N6,bodyX,18) (N6,bodyY,19)
            // (N7,tip0x,20) (N7,tip0y,21) (N7,tip1x,22) (N7,tip1y,23) (N7,tip2x,24) (N7,tip2y,25) (N7,tip3x,26) (N7,tip3y,27)

            // Calculate relevant stiffness submatrices from elements
            IMatrix elem1Stiffness = ((XContinuumElement2D)model.Elements[1].ElementType).JoinStiffnessesStandardFirst();
            IMatrix elem2Stiffness = ((XContinuumElement2D)model.Elements[2].ElementType).JoinStiffnessesStandardFirst();

            int[]   elem2Node6Dofs = { 4, 5, 18, 19 };
            IMatrix node6Stiffness = elem2Stiffness.GetSubmatrix(elem2Node6Dofs, elem2Node6Dofs);

            int[]   elem1Node7Dofs      = { 4, 5, 24, 25, 26, 27, 28, 29, 30, 31 };
            IMatrix node7Elem1Stiffness = elem1Stiffness.GetSubmatrix(elem1Node7Dofs, elem1Node7Dofs);

            int[]   elem2Node7Dofs      = { 6, 7, 20, 21, 22, 23, 24, 25, 26, 27 };
            IMatrix node7Elem2Stiffness = elem2Stiffness.GetSubmatrix(elem2Node7Dofs, elem2Node7Dofs);

            // Check matrices
            double equalityTolerance    = 1E-13;
            Func <double, double> round = x => 1E6 * Math.Round(x * 1E-6, 3);

            Assert.True(node6StiffnessExpected.Equals(node6Stiffness.DoToAllEntries(round), equalityTolerance));
            Assert.True(node7Elem1StiffnessExpected.Equals(node7Elem1Stiffness.DoToAllEntries(round), equalityTolerance));
            Assert.True(node7Elem2StiffnessExpected.Equals(node7Elem2Stiffness.DoToAllEntries(round), equalityTolerance));
            Assert.True(node7GlobalStiffnessExpected.Equals(node7GlobalStiffness.DoToAllEntries(round), equalityTolerance));
        }