/usr/share/doc/freefem++/examples/examples++-tutorial/sparse-cmatrix.edp is in freefem++-doc 3.19.1-1.
This file is owned by root:root, with mode 0o644.
The actual contents of the file can be viewed below.
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 | // sparse matrix test ---
// example of the new matrix feature in version 1.40
// -------------------------------------------------
mesh TH = square(3,4);
mesh th = square(2,3);
mesh Th = square(4,4);
complex ccc;
ccc= 1;
cout << ccc << endl;
fespace VH(TH,P1);
fespace Vh(th,P1);
fespace Wh(Th,P1);
matrix RB= interpolate(VH,Vh); // build interpolation matrix Vh->Vh
matrix RBB= interpolate(Wh,Vh); // build interpolation matrix
matrix<complex> B=RB;
B = B*(1+2i);
matrix<complex> BB=RBB;
varf vA(u,v) = int2d(Th)(dx(u)*dx(v)+dy(u)*dy(v))+ int1d(Th)(u*v);
matrix<complex> A=vA(Wh,Wh);
Vh<complex> ml=0;
cout << " ml " << ml[] << endl;
varf vML(u,v) = int2d(th)(1.*v);
ml[]=vML(0,Vh); // build the P1 mass lump of P1
cout << ml[] << endl;
matrix<complex> ML(ml[]); // matrix diagonal
cout << "ML="<<ML << endl;
cout << "B="<<B << endl;
matrix<complex> BML=B*ML; // a faire
matrix<complex> tB=B'; //'; transpose and conjugate
cout << "tB=" << tB << endl;
matrix<complex> MLtB=ML'*B'; //
cout << "BML="<<BML << endl;
cout << "MLtB=" << MLtB << endl;
// WARNING if UMFPACK is not install
// the UMFPACK solver is replace by LU
// but LU need skyline matrix
int typesolver=UMFPACK;
if(HaveUMFPACK) typesolver=GMRES;
set(A,solver=typesolver); // set a solver
VH<complex> uH=0;
Vh<complex> uh=x+y+1i*(x-y);
uH[]= B*uh[];
Vh uHr = imag(uH);
plot(uHr,wait=1);
matrix<complex> BtA = BB'*A;//';
matrix<complex> BtAB = BtA*BB;
if(HaveUMFPACK)
set(BtAB,solver=UMFPACK);
else
set(BtAB,solver=GMRES);
Vh<complex> ff=1+1i;
Vh<complex> xx;
Vh xxr;
cout << " ------ " << endl;
xx[]=BtAB^-1*ff[];
cout << " ------ " << endl;
xx[]=BtAB^-1*ff[];
cout << " ------ " << endl;
xxr=imag(xx);
plot(xxr, wait=1);
{
int N=10;
complex [int,int] A(N,N);
complex [int] a(N),b(N),bb(N);
A =0;
for (int i=0;i<N;i++)
{
A(i,i)=1.+i;
if(i+1 < N) A(i,i+1)=-i-1i*i;
a[i]=i*(1.+2i);
}
b=A*a;
cout << " b =" << b << endl ;
cout << " a =" << a << endl ;
cout << " b'*b (hermissian product) = " << b'*b << endl; //';
cout << " a'*a = " << a'*a << endl;//';;
assert( abs(imag(b'*b)) <1e-5);//')));
cout << "xxxx\n";
matrix<complex> sparseA=A;
cout << sparseA << endl;
sparseA = 2*sparseA+sparseA;
sparseA = 4*sparseA+sparseA*(5+1i); // * 27
matrix<complex> sparseB=sparseA;//+sparseA+sparseA; ;
cout << sparseA << endl;
cout << sparseB << endl; // *81
cout << "sparseB = " << sparseB(0,0) << endl;
// ajoute version 2.0-2
sparseA=A;
verbosity=4;
if(HaveUMFPACK)
set(sparseA,solver=UMFPACK,tolpivot=1e-10,tolpivotsym=1e-9);
else
set(sparseA,solver=GMRES);
bb=sparseA^-1*a;
verbosity=1;
b = sparseA*bb;
b -= a;
cout << " nb coef sparseA " << sparseA.nbcoef << endl;
cout << " ||b.||_1 " << b.l1 << endl;
cout << " ||b.||_2 " << b.l2 << endl;
cout << " ||b.||_infty " << b.linfty << endl;
assert(b.l1 < 1e-10);
}
{// version 3.8
mesh Th=square(2,2);
fespace Xh(Th,P1);
varf vv(u,v)= int2d(Th)( ((2+1i)*u*v)')+int2d(Th)((3+2i)*v);//');
varf vr(u,v)= int2d(Th)( u*v);//');
matrix<complex> A=vv(Xh,Xh);
matrix Ar=vr(Xh,Xh);
complex[int] vc=vv(0,Xh);
real[int] vrr=vc.re,vii=vc.im;
vrr=vc.re;
vii=vc.im;
cout << "vc[0] = " <<vc[0] << " = " << vc.re[0] << " +i " << vc.im[0] <<endl;
cout << [ 1i, 1i]'*[ 1i, 1i] <<endl;//';
// real part un complex par of matrix .
Ar = A.re;
cout <<" A(0,0) = " << A(0,0) << " ";
cout << " A.re(0.0) = " << Ar(0,0) << " " ;
Ar = A.im;
cout << " A.im(0.0) = " << Ar(0,0) << endl ;
macro Grada(u) [ phia*dx(u) + phiax*u ,dy(u) ]// ...
func phia = exp(-2i*pi*x);
func phiax = -2i*pi*exp(-2i*pi*x);
varf va(u,v)= int2d(Th)( Grada(v)'*Grada(u) ) ;//');
A = va(Xh,Xh);
matrix<complex> At=A';
cout << A(1,2)' << " == " << At(2,1) << endl;
A = A+ (-1)*At;
cout << A << endl;
// copy and initialisation of complex matric with real matrix.
A = Ar;
matrix<complex> Ac=Ar;
}
|