summaryrefslogtreecommitdiff
path: root/src/FEmTool/FEmTool_LinearTension.cxx
blob: b092b552b95b7125e8a28488362f4c2edbb831f2 (plain)
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
// File:	FEmTool_LinearTension.cxx
// Created:	Fri Nov  6 10:22:02 1998
// Author:	Igor FEOKTISTOV
//		<ifv@paradox.nnov.matra-dtv.fr>


#include <FEmTool_LinearTension.ixx>
#include <PLib.hxx>
#include <TColStd_HArray2OfInteger.hxx>
#include <TColStd_HArray2OfReal.hxx>
#include <PLib_JacobiPolynomial.hxx>
#include <PLib_HermitJacobi.hxx>
#include <FEmTool_ElementsOfRefMatrix.hxx>
#include <math_IntegerVector.hxx>
#include <math_Vector.hxx>
#include <math_GaussSetIntegration.hxx>
#include <math.hxx>
#include <Standard_ConstructionError.hxx>

FEmTool_LinearTension::FEmTool_LinearTension(const Standard_Integer WorkDegree,
					     const GeomAbs_Shape ConstraintOrder):
       RefMatrix(0,WorkDegree,0,WorkDegree)
       
{
  static Standard_Integer Order = -333, WDeg = 14;
  static math_Vector MatrixElemts(0, ((WDeg+2)*(WDeg+1))/2 -1 );

  myOrder = PLib::NivConstr(ConstraintOrder);

  if (myOrder != Order) {
    //Calculating RefMatrix
    if (WorkDegree > WDeg) Standard_ConstructionError::Raise("Degree too high");
    Order = myOrder;
    Standard_Integer DerOrder = 1;
    Handle(PLib_HermitJacobi) theBase = new PLib_HermitJacobi(WDeg, ConstraintOrder);
    FEmTool_ElementsOfRefMatrix Elem = FEmTool_ElementsOfRefMatrix(theBase, DerOrder);
    
    Standard_Integer maxDegree = WDeg+1;
    math_IntegerVector Order(1,1,Min(4*(maxDegree/2+1),math::GaussPointsMax()));
    math_Vector Lower(1,1,-1.), Upper(1,1,1.); 
    
    math_GaussSetIntegration anInt(Elem, Lower, Upper, Order);
    MatrixElemts = anInt.Value();
  }

  Standard_Integer i, j, ii, jj;
  for(ii = i = 0; i <= WorkDegree; i++) {
    RefMatrix(i, i) =  MatrixElemts(ii);
    for(j = i+1, jj = ii+1; j <= WorkDegree; j++, jj++) {
      RefMatrix(j, i) = RefMatrix(i, j) =  MatrixElemts(jj);
    }
    ii += WDeg+1-i;
  }
}

Handle(TColStd_HArray2OfInteger) FEmTool_LinearTension::DependenceTable() const
{
  if(myCoeff.IsNull()) Standard_DomainError::Raise("FEmTool_LinearTension::DependenceTable");

  Handle(TColStd_HArray2OfInteger) DepTab = 
    new TColStd_HArray2OfInteger(myCoeff->LowerCol(), myCoeff->UpperCol(),
				 myCoeff->LowerCol(), myCoeff->UpperCol(),0);
  Standard_Integer i;
  for(i=1; i<=myCoeff->RowLength(); i++) DepTab->SetValue(i,i,1);
  
  return DepTab;
}

Standard_Real FEmTool_LinearTension::Value() 
{
  Standard_Integer deg = Min(myCoeff->ColLength() - 1, RefMatrix.UpperRow()), 
                   i, j, j0 = myCoeff->LowerRow(), degH = Min(2*myOrder+1, deg),
                   NbDim = myCoeff->RowLength(), dim;

  TColStd_Array2OfReal NewCoeff( 1, NbDim, 0, deg);

  Standard_Real coeff = (myLast - myFirst)/2., cteh3 = 2./coeff, 
                mfact, Jline;

  Standard_Integer k1;


  Standard_Real J = 0.;

  for(i = 0; i <= degH; i++) {
    k1 = (i <= myOrder)? i : i - myOrder - 1;
    mfact = Pow(coeff,k1);
    for(dim = 1; dim <= NbDim; dim++) 
      NewCoeff(dim, i) = myCoeff->Value(j0 + i, dim) * mfact;
  }

  for(i = degH + 1; i <= deg; i++) {
    for(dim = 1; dim <= NbDim; dim++) 
      NewCoeff(dim, i) = myCoeff->Value(j0 + i, dim);
  }

  for(dim = 1; dim <= NbDim; dim++) {
  
    for(i = 0; i <= deg; i++) {

      Jline = 0.5 * RefMatrix(i, i) * NewCoeff(dim, i);

      for(j = 0; j < i; j++) 
	Jline += RefMatrix(i, j) * NewCoeff(dim, j);
      
      J += Jline * NewCoeff(dim, i);

    }    

  }

  return cteh3*J;


}

void FEmTool_LinearTension::Hessian(const Standard_Integer Dimension1,
				    const Standard_Integer Dimension2, math_Matrix& H) 
{
 
  Handle(TColStd_HArray2OfInteger) DepTab = DependenceTable();

  if(Dimension1 < DepTab->LowerRow() || Dimension1 > DepTab->UpperRow() || 
     Dimension2 < DepTab->LowerCol() || Dimension2 > DepTab->UpperCol()) 
    Standard_OutOfRange::Raise("FEmTool_LinearTension::Hessian");

  if(DepTab->Value(Dimension1,Dimension2) == 0) 
    Standard_DomainError::Raise("FEmTool_LinearTension::Hessian");

  Standard_Integer deg = Min(RefMatrix.UpperRow(), H.RowNumber() - 1), degH = Min(2*myOrder+1, deg);

  Standard_Real coeff = (myLast - myFirst)/2., cteh3 = 2./coeff, mfact;
  Standard_Integer k1, k2, i, j, i0 = H.LowerRow(), j0 = H.LowerCol(), i1, j1;

  H.Init(0.);

  i1 = i0;
  for(i = 0; i <= degH; i++) {
    k1 = (i <= myOrder)? i : i - myOrder - 1;
    mfact = Pow(coeff,k1)*cteh3;
    // Hermite*Hermite part of matrix
    j1 = j0 + i;
    for(j = i; j <= degH; j++) {
      k2 = (j <= myOrder)? j : j - myOrder - 1;
      H(i1, j1) = mfact*Pow(coeff, k2)*RefMatrix(i, j);
      if (i != j) H(j1, i1) = H(i1, j1);
      j1++;
    }
    // Hermite*Jacobi part of matrix
    j1 = j0 + degH + 1;
    for(j = degH + 1; j <= deg; j++) {
      H(i1, j1) = mfact*RefMatrix(i, j);
      H(j1, i1) = H(i1, j1);
      j1++;
    }
    i1++;
  }


  // Jacoby*Jacobi part of matrix
  i1 = i0 + degH + 1;
  for(i = degH+1; i <= deg; i++) {
    j1 = j0 + i;
    for(j = i; j <= deg; j++) {
      H(i1, j1) = cteh3*RefMatrix(i, j);
      if (i != j) H(j1, i1) = H(i1, j1);
      j1++;
    }
    i1++;
  }

}

 void FEmTool_LinearTension::Gradient(const Standard_Integer Dimension, math_Vector& G) 
{
  if(Dimension < myCoeff->LowerCol() || Dimension > myCoeff->UpperCol()) 
    Standard_OutOfRange::Raise("FEmTool_LinearTension::Gradient");

  Standard_Integer deg = Min(G.Length() - 1, myCoeff->ColLength() - 1);

  math_Vector X(0,deg);
  Standard_Integer i, i1 = myCoeff->LowerRow();
  for(i = 0; i <= deg; i++) X(i) = myCoeff->Value(i1+i, Dimension);

  math_Matrix H(0,deg,0,deg);
  Hessian(Dimension, Dimension, H);
  
  G.Multiply(H, X);


}