Transferred copyright to the OpenFOAM Foundation
[OpenFOAM-2.0.x.git] / applications / solvers / electromagnetics / magneticFoam / createFields.H
blob22d40089213178d452b55a42c2b3ecd503afa61d
1     Info<< "Reading field psi\n" << endl;
2     volScalarField psi
3     (
4         IOobject
5         (
6             "psi",
7             runTime.timeName(),
8             mesh,
9             IOobject::MUST_READ,
10             IOobject::AUTO_WRITE
11         ),
12         mesh
13     );
15     Info<< "Reading transportProperties\n" << endl;
17     IOdictionary transportProperties
18     (
19         IOobject
20         (
21             "transportProperties",
22             runTime.constant(),
23             mesh,
24             IOobject::MUST_READ_IF_MODIFIED,
25             IOobject::NO_WRITE
26         )
27     );
29     List<magnet> magnets(transportProperties.lookup("magnets"));
31     surfaceScalarField murf
32     (
33         IOobject
34         (
35             "murf",
36             runTime.timeName(),
37             mesh
38         ),
39         mesh,
40         1
41     );
43     surfaceScalarField Mrf
44     (
45         IOobject
46         (
47             "Mrf",
48             runTime.timeName(),
49             mesh
50         ),
51         mesh,
52         dimensionedScalar("Mr", dimensionSet(0, 1, 0, 0, 0, 1, 0), 0)
53     );
55     forAll(magnets, i)
56     {
57         label magnetZonei = mesh.faceZones().findZoneID(magnets[i].name());
59         if (magnetZonei == -1)
60         {
61             FatalIOErrorIn(args.executable().c_str(), transportProperties)
62                 << "Cannot find faceZone for magnet " << magnets[i].name()
63                 << exit(FatalIOError);
64         }
66         const labelList& faces = mesh.faceZones()[magnetZonei];
68         const scalar muri = magnets[i].mur();
69         const scalar Mri = magnets[i].Mr().value();
70         const vector& orientationi = magnets[i].orientation();
72         const surfaceVectorField& Sf = mesh.Sf();
74         forAll(faces, i)
75         {
76             label facei = faces[i];
77             murf[facei] = muri;
78             Mrf[facei] = Mri*(orientationi & Sf[facei]);
79         }
80     }