00001
00050 #include "apbscfg.h"
00051
00052 #ifdef HAVE_MC_H
00053
00054 #include "apbs/vfetk.h"
00055
00056
00057
00058
00059
00060
00061
00062
00063
00064
00065
00066
00067
00068
00069
00070
00071 VPRIVATE double Vfetk_qfEnergyAtom(
00072 Vfetk *thee,
00073 int iatom,
00074 int color,
00075 double *sol
00076 );
00077
00078
00079
00080
00081
00082
00083 VPRIVATE Vfetk_LocalVar var;
00084
00085
00086
00087
00088
00089
00090 VPRIVATE char *diriCubeString =
00091 "mcsf_begin=1;\n\
00092 \n\
00093 dim=3;\n\
00094 dimii=3;\n\
00095 vertices=8;\n\
00096 simplices=6;\n\
00097 \n\
00098 vert=[\n\
00099 0 0 -0.5 -0.5 -0.5\n\
00100 1 0 0.5 -0.5 -0.5\n\
00101 2 0 -0.5 0.5 -0.5\n\
00102 3 0 0.5 0.5 -0.5\n\
00103 4 0 -0.5 -0.5 0.5\n\
00104 5 0 0.5 -0.5 0.5\n\
00105 6 0 -0.5 0.5 0.5\n\
00106 7 0 0.5 0.5 0.5\n\
00107 ];\n\
00108 \n\
00109 simp=[\n\
00110 0 0 0 0 1 0 1 0 5 1 2\n\
00111 1 0 0 0 1 1 0 0 5 2 4\n\
00112 2 0 0 0 1 0 1 1 5 3 2\n\
00113 3 0 0 0 1 0 1 3 5 7 2\n\
00114 4 0 0 1 1 0 0 2 5 7 6\n\
00115 5 0 0 1 1 0 0 2 5 6 4\n\
00116 ];\n\
00117 \n\
00118 mcsf_end=1;\n\
00119 \n\
00120 ";
00121
00122
00123
00124
00125
00126
00127 VPRIVATE char *neumCubeString =
00128 "mcsf_begin=1;\n\
00129 \n\
00130 dim=3;\n\
00131 dimii=3;\n\
00132 vertices=8;\n\
00133 simplices=6;\n\
00134 \n\
00135 vert=[\n\
00136 0 0 -0.5 -0.5 -0.5\n\
00137 1 0 0.5 -0.5 -0.5\n\
00138 2 0 -0.5 0.5 -0.5\n\
00139 3 0 0.5 0.5 -0.5\n\
00140 4 0 -0.5 -0.5 0.5\n\
00141 5 0 0.5 -0.5 0.5\n\
00142 6 0 -0.5 0.5 0.5\n\
00143 7 0 0.5 0.5 0.5\n\
00144 ];\n\
00145 \n\
00146 simp=[\n\
00147 0 0 0 0 2 0 2 0 5 1 2\n\
00148 1 0 0 0 2 2 0 0 5 2 4\n\
00149 2 0 0 0 2 0 2 1 5 3 2\n\
00150 3 0 0 0 2 0 2 3 5 7 2\n\
00151 4 0 0 2 2 0 0 2 5 7 6\n\
00152 5 0 0 2 2 0 0 2 5 6 4\n\
00153 ];\n\
00154 \n\
00155 mcsf_end=1;\n\
00156 \n\
00157 ";
00158
00159
00160
00161
00162
00163
00164
00165
00166
00167 VPRIVATE double diel();
00168
00169
00170
00171
00172
00173
00174
00175
00176
00177 VPRIVATE double ionacc();
00178
00179
00180
00181
00182
00183
00184
00185
00186
00187
00188
00189
00190
00191 VPRIVATE double smooth(
00192 int nverts,
00193 double dist[VAPBS_NVS],
00194 double coeff[VAPBS_NVS],
00195 int meth
00196 );
00197
00198
00199
00200
00201
00202
00203
00204
00205
00206
00207
00208
00209
00210 VPRIVATE double debye_U(
00211 Vpbe *pbe,
00212 int d,
00213 double x[]
00214 );
00215
00216
00217
00218
00219
00220
00221
00222
00223
00224
00225
00226
00227 VPRIVATE double debye_Udiff(
00228 Vpbe *pbe,
00229 int d,
00230 double x[]
00231 );
00232
00233
00234
00235
00236
00237
00238
00239
00240
00241
00242
00243
00244
00245
00246
00247 VPRIVATE void coulomb(
00248 Vpbe *pbe,
00249 int d,
00250 double x[],
00251 double eps,
00252 double *U,
00253 double dU[],
00254 double *d2U
00255 );
00256
00257
00258
00259
00260
00261
00262
00263
00264
00265
00266
00267 VPRIVATE void init_2DP1(
00268 int dimIS[],
00269 int *ndof,
00270 int dof[],
00271 double c[][VMAXP],
00272 double cx[][VMAXP],
00273 double cy[][VMAXP],
00274 double cz[][VMAXP]
00275 );
00276
00277
00278
00279
00280
00281
00282
00283
00284
00285
00286
00287
00288
00289 VPRIVATE void init_3DP1(
00290 int dimIS[],
00291 int *ndof,
00292 int dof[],
00293 double c[][VMAXP],
00294 double cx[][VMAXP],
00295 double cy[][VMAXP],
00296 double cz[][VMAXP]
00297 );
00298
00299
00300
00301
00302
00303
00304
00305
00306
00307
00308
00309
00310
00311
00312
00313 VPRIVATE void setCoef(
00314 int numP,
00315 double c[][VMAXP],
00316 double cx[][VMAXP],
00317 double cy[][VMAXP],
00318 double cz[][VMAXP],
00319 int ic[][VMAXP],
00320 int icx[][VMAXP],
00321 int icy[][VMAXP],
00322 int icz[][VMAXP]
00323 );
00324
00325
00326
00327
00328
00329
00330
00331
00332
00333
00334
00335
00336
00337
00338
00339
00340
00341
00342
00343
00344
00345
00346
00347
00348
00349
00350
00351 VPRIVATE void polyEval(
00352 int numP,
00353 double p[],
00354 double c[][VMAXP],
00355 double xv[]
00356 );
00357
00358
00359
00360
00361
00362
00363 VPRIVATE int dim_2DP1 = 3;
00364
00365
00366
00367
00368
00369
00370
00371
00372
00373
00374
00375
00376
00377
00378
00379
00380
00381
00382
00383 VPRIVATE int lgr_2DP1[3][VMAXP] = {
00384
00385
00386
00387
00388 { 2, -2, -2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 },
00389 { 0, 2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 },
00390 { 0, 0, 2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 }
00391 };
00392 VPRIVATE int lgr_2DP1x[3][VMAXP] = {
00393
00394
00395 { -2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 },
00396 { 2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 },
00397 { 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 }
00398 };
00399 VPRIVATE int lgr_2DP1y[3][VMAXP] = {
00400
00401
00402 { -2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 },
00403 { 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 },
00404 { 2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 }
00405 };
00406 VPRIVATE int lgr_2DP1z[3][VMAXP] = {
00407
00408
00409 { 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 },
00410 { 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 },
00411 { 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 }
00412 };
00413
00414
00415
00416
00417
00418
00419
00420
00421
00422
00423
00424
00425
00426
00427
00428
00429
00430
00431
00432
00433
00434 VPRIVATE int dim_3DP1 = VAPBS_NVS;
00435 VPRIVATE int lgr_3DP1[VAPBS_NVS][VMAXP] = {
00436
00437
00438 { 2, -2, -2, -2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 },
00439 { 0, 2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 },
00440 { 0, 0, 2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 },
00441 { 0, 0, 0, 2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 }
00442 };
00443 VPRIVATE int lgr_3DP1x[VAPBS_NVS][VMAXP] = {
00444
00445
00446 { -2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 },
00447 { 2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 },
00448 { 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 },
00449 { 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 }
00450 };
00451 VPRIVATE int lgr_3DP1y[VAPBS_NVS][VMAXP] = {
00452
00453
00454 { -2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 },
00455 { 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 },
00456 { 2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 },
00457 { 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 }
00458 };
00459 VPRIVATE int lgr_3DP1z[VAPBS_NVS][VMAXP] = {
00460
00461
00462 { -2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 },
00463 { 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 },
00464 { 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 },
00465 { 2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 },
00466 };
00467
00468
00469
00470
00471
00472
00473 VPRIVATE const int P_DEG=1;
00474
00475
00476
00477
00478
00479 VPRIVATE int numP;
00480 VPRIVATE double c[VMAXP][VMAXP];
00481 VPRIVATE double cx[VMAXP][VMAXP];
00482 VPRIVATE double cy[VMAXP][VMAXP];
00483 VPRIVATE double cz[VMAXP][VMAXP];
00484
00485 #if !defined(VINLINE_VFETK)
00486
00487 VPUBLIC Gem* Vfetk_getGem(Vfetk *thee) {
00488
00489 VASSERT(thee != VNULL);
00490 return thee->gm;
00491
00492 }
00493
00494 VPUBLIC AM* Vfetk_getAM(Vfetk *thee) {
00495
00496 VASSERT(thee != VNULL);
00497 return thee->am;
00498 }
00499
00500 VPUBLIC Vpbe* Vfetk_getVpbe(Vfetk *thee) {
00501
00502 VASSERT(thee != VNULL);
00503 return thee->pbe;
00504
00505 }
00506
00507 VPUBLIC Vcsm* Vfetk_getVcsm(Vfetk *thee) {
00508
00509 VASSERT(thee != VNULL);
00510 return thee->csm;
00511
00512 }
00513
00514 VPUBLIC int Vfetk_getAtomColor(Vfetk *thee, int iatom) {
00515
00516 int natoms;
00517
00518 VASSERT(thee != VNULL);
00519
00520 natoms = Valist_getNumberAtoms(Vpbe_getValist(thee->pbe));
00521 VASSERT(iatom < natoms);
00522
00523 return Vatom_getPartID(Valist_getAtom(Vpbe_getValist(thee->pbe), iatom));
00524 }
00525 #endif
00526
00527 VPUBLIC Vfetk* Vfetk_ctor(Vpbe *pbe, Vhal_PBEType type) {
00528
00529
00530 Vfetk *thee = VNULL;
00531 thee = Vmem_malloc(VNULL, 1, sizeof(Vfetk) );
00532 VASSERT(thee != VNULL);
00533 VASSERT(Vfetk_ctor2(thee, pbe, type));
00534
00535 return thee;
00536 }
00537
00538 VPUBLIC int Vfetk_ctor2(Vfetk *thee, Vpbe *pbe, Vhal_PBEType type) {
00539
00540 int i;
00541 double center[VAPBS_DIM];
00542
00543
00544 VASSERT(pbe != VNULL);
00545 thee->pbe = pbe;
00546 VASSERT(pbe->alist != VNULL);
00547 VASSERT(pbe->acc != VNULL);
00548
00549
00550 thee->type = type;
00551
00552
00553 thee->vmem = Vmem_ctor("APBS::VFETK");
00554
00555
00556 Vnm_print(0, "Vfetk_ctor2: Constructing PDE...\n");
00557 thee->pde = Vfetk_PDE_ctor(thee);
00558 Vnm_print(0, "Vfetk_ctor2: Constructing Gem...\n");
00559 thee->gm = Gem_ctor(thee->vmem, thee->pde);
00560 Vnm_print(0, "Vfetk_ctor2: Constructing Aprx...\n");
00561 thee->aprx = Aprx_ctor(thee->vmem, thee->gm, thee->pde);
00562 Vnm_print(0, "Vfetk_ctor2: Constructing Aprx...\n");
00563 thee->am = AM_ctor(thee->vmem, thee->aprx);
00564
00565
00566 thee->level = 0;
00567
00568
00569 thee->lkey = VLT_MG;
00570 thee->lmax = 1000000;
00571 thee->ltol = 1e-5;
00572 thee->lprec = VPT_MG;
00573 thee->nkey = VNT_NEW;
00574 thee->nmax = 1000000;
00575 thee->ntol = 1e-5;
00576 thee->gues = VGT_ZERO;
00577 thee->pjac = -1;
00578
00579
00580 var.fetk = thee;
00581 var.initGreen = 0;
00582
00583
00584 Gem_setExternalUpdateFunction(thee->gm, Vfetk_externalUpdateFunction);
00585
00586
00587 var.zkappa2 = Vpbe_getZkappa2(var.fetk->pbe);
00588 var.ionstr = Vpbe_getBulkIonicStrength(var.fetk->pbe);
00589 if (var.ionstr > 0.0) var.zks2 = 0.5*var.zkappa2/var.ionstr;
00590 else var.zks2 = 0.0;
00591 Vpbe_getIons(var.fetk->pbe, &(var.nion), var.ionConc, var.ionRadii,
00592 var.ionQ);
00593 for (i=0; i<var.nion; i++) {
00594 var.ionConc[i] = var.zks2 * var.ionConc[i] * var.ionQ[i];
00595 }
00596
00597
00598 thee->pbeparm = VNULL;
00599 thee->feparm = VNULL;
00600 thee->csm = VNULL;
00601
00602 return 1;
00603 }
00604
00605 VPUBLIC void Vfetk_setParameters(Vfetk *thee, PBEparm *pbeparm,
00606 FEMparm *feparm) {
00607
00608 VASSERT(thee != VNULL);
00609 thee->feparm = feparm;
00610 thee->pbeparm = pbeparm;
00611 }
00612
00613 VPUBLIC void Vfetk_dtor(Vfetk **thee) {
00614 if ((*thee) != VNULL) {
00615 Vfetk_dtor2(*thee);
00616
00617 (*thee) = VNULL;
00618 }
00619 }
00620
00621 VPUBLIC void Vfetk_dtor2(Vfetk *thee) {
00622 Vcsm_dtor(&(thee->csm));
00623 AM_dtor(&(thee->am));
00624 Aprx_dtor(&(thee->aprx));
00625 Vfetk_PDE_dtor(&(thee->pde));
00626 Vmem_dtor(&(thee->vmem));
00627 }
00628
00629 VPUBLIC double* Vfetk_getSolution(Vfetk *thee, int *length) {
00630
00631 int i;
00632 double *solution;
00633 double *theAnswer;
00634 AM *am;
00635
00636 VASSERT(thee != VNULL);
00637
00638
00639 am = thee->am;
00640
00641 Bvec_copy(am->w0, am->u);
00642
00643 Bvec_axpy(am->w0, am->ud, 1.);
00644
00645 solution = Bvec_addr(am->w0);
00646
00647 *length = Bvec_numRT(am->w0);
00648
00649
00650 VASSERT(1 == Bvec_numB(am->w0));
00651
00652 theAnswer = VNULL;
00653 theAnswer = Vmem_malloc(VNULL, *length, sizeof(double));
00654 VASSERT(theAnswer != VNULL);
00655 for (i=0; i<(*length); i++) theAnswer[i] = solution[i];
00656
00657 return theAnswer;
00658 }
00659
00660 VPUBLIC double Vfetk_energy(Vfetk *thee, int color, int nonlin) {
00661
00662 double totEnergy = 0.0;
00663 double qfEnergy = 0.0;
00664 double dqmEnergy = 0.0;
00665
00666 VASSERT(thee != VNULL);
00667
00668 if (nonlin && (Vpbe_getBulkIonicStrength(thee->pbe) > 0.)) {
00669 Vnm_print(0, "Vfetk_energy: calculating full PBE energy\n");
00670 Vnm_print(0, "Vfetk_energy: bulk ionic strength = %g M\n",
00671 Vpbe_getBulkIonicStrength(thee->pbe));
00672 dqmEnergy = Vfetk_dqmEnergy(thee, color);
00673 Vnm_print(0, "Vfetk_energy: dqmEnergy = %g kT\n", dqmEnergy);
00674 qfEnergy = Vfetk_qfEnergy(thee, color);
00675 Vnm_print(0, "Vfetk_energy: qfEnergy = %g kT\n", qfEnergy);
00676
00677 totEnergy = qfEnergy - dqmEnergy;
00678 } else {
00679 Vnm_print(0, "Vfetk_energy: calculating only q-phi energy\n");
00680 dqmEnergy = Vfetk_dqmEnergy(thee, color);
00681 Vnm_print(0, "Vfetk_energy: dqmEnergy = %g kT (NOT USED)\n", dqmEnergy);
00682 qfEnergy = Vfetk_qfEnergy(thee, color);
00683 Vnm_print(0, "Vfetk_energy: qfEnergy = %g kT\n", qfEnergy);
00684 totEnergy = 0.5*qfEnergy;
00685 }
00686
00687 return totEnergy;
00688
00689 }
00690
00691
00692 VPUBLIC double Vfetk_qfEnergy(Vfetk *thee, int color) {
00693
00694 double *sol; int nsol;
00695 int iatom, natoms;
00696 double energy = 0.0;
00697
00698 AM *am;
00699
00700 VASSERT(thee != VNULL);
00701 am = thee->am;
00702
00703
00704 sol= VNULL;
00705 sol = Vfetk_getSolution(thee, &nsol);
00706 VASSERT(sol != VNULL);
00707
00708
00709
00710 if (nsol != Gem_numVV(thee->gm)) {
00711 Vnm_print(2, "Vfetk_qfEnergy: Number of unknowns in solution does not match\n");
00712 Vnm_print(2, "Vfetk_qfEnergy: number of vertices in mesh!!! Bailing out!\n");
00713 VASSERT(0);
00714 }
00715
00716
00717 natoms = Valist_getNumberAtoms(thee->pbe->alist);
00718 for (iatom=0; iatom<natoms; iatom++) {
00719
00720 energy = energy + Vfetk_qfEnergyAtom(thee, iatom, color, sol);
00721
00722 }
00723
00724
00725 Vmem_free(VNULL, nsol, sizeof(double), (void **)&sol);
00726
00727
00728 return energy;
00729 }
00730
00731 VPRIVATE double Vfetk_qfEnergyAtom(
00732 Vfetk *thee,
00733 int iatom,
00734 int color,
00735 double *sol) {
00736
00737 Vatom *atom;
00738 double charge;
00739 double phi[VAPBS_NVS], phix[VAPBS_NVS][3], *position;
00740 double uval;
00741 double energy = 0.0;
00742 int isimp, nsimps;
00743 SS *simp;
00744 int icolor, ivert, usingColor;
00745
00746
00747
00748 atom = Valist_getAtom(thee->pbe->alist, iatom);
00749 icolor = Vfetk_getAtomColor(thee, iatom);
00750 charge = Vatom_getCharge(atom);
00751 position = Vatom_getPosition(atom);
00752
00753
00754 usingColor = (color >= 0);
00755
00756 if (usingColor && (icolor<0)) {
00757 Vnm_print(2, "Vfetk_qfEnergy: Atom colors not set!\n");
00758 VASSERT(0);
00759 }
00760
00761
00762 if ((icolor==color) || (!usingColor)) {
00763
00764 nsimps = Vcsm_getNumberSimplices(thee->csm, iatom);
00765
00766
00767
00768
00769 for (isimp=0; isimp<nsimps; isimp++) {
00770
00771 simp = Vcsm_getSimplex(thee->csm, isimp, iatom);
00772
00773
00774
00775 if ((SS_chart(simp)==color)||(color<0)) {
00776
00777
00778 Gem_pointInSimplexVal(thee->gm, simp, position, phi, phix);
00779 for (ivert=0; ivert<SS_dimVV(simp); ivert++) {
00780 uval = sol[VV_id(SS_vertex(simp,ivert))];
00781 energy += (charge*phi[ivert]*uval);
00782 }
00783
00784
00785 break;
00786 }
00787 }
00788 }
00789
00790 return energy;
00791 }
00792
00793
00794 VPUBLIC double Vfetk_dqmEnergy(Vfetk *thee, int color) {
00795
00796 return AM_evalJ(thee->am);
00797
00798 }
00799
00800 VPUBLIC void Vfetk_setAtomColors(Vfetk *thee) {
00801
00802 #define VMAXLOCALCOLORSDONTREUSETHISVARIABLE 1024
00803 SS *simp;
00804 Vatom *atom;
00805 int i, natoms;
00806
00807 VASSERT(thee != VNULL);
00808
00809 natoms = Valist_getNumberAtoms(thee->pbe->alist);
00810 for (i=0; i<natoms; i++) {
00811 atom = Valist_getAtom(thee->pbe->alist, i);
00812 simp = Vcsm_getSimplex(thee->csm, 0, i);
00813 Vatom_setPartID(atom, SS_chart(simp));
00814 }
00815
00816 }
00817
00818 VPUBLIC unsigned long int Vfetk_memChk(Vfetk *thee) {
00819
00820 int memUse = 0;
00821
00822 if (thee == VNULL) return 0;
00823
00824 memUse = memUse + sizeof(Vfetk);
00825 memUse = memUse + Vcsm_memChk(thee->csm);
00826
00827 return memUse;
00828 }
00829
00830 VPUBLIC Vrc_Codes Vfetk_genCube(
00831 Vfetk *thee,
00832 double center[3],
00833 double length[3],
00834 Vfetk_MeshLoad meshType) {
00835 AM *am = VNULL;
00836 Gem *gm = VNULL;
00837
00838 int skey = 0;
00839 char *key = "r";
00840 char *iodev = "BUFF";
00841 char *iofmt = "ASC";
00842 char *iohost = "localhost";
00843 char *iofile = "0";
00844 Vio *sock = VNULL;
00845 char buf[VMAX_BUFSIZE];
00846 int bufsize = 0;
00847 VV *vx = VNULL;
00848 int i, j;
00849 double x;
00850
00851 VASSERT(thee != VNULL);
00852 am = thee->am;
00853 VASSERT(am != VNULL);
00854 gm = thee->gm;
00855 VASSERT(gm != VNULL);
00856
00857
00858
00859 switch (meshType) {
00860 case VML_DIRICUBE:
00861 bufsize = strlen(diriCubeString);
00862 VASSERT( bufsize <= VMAX_BUFSIZE );
00863 strncpy(buf, diriCubeString, VMAX_BUFSIZE);
00864 break;
00865 case VML_NEUMCUBE:
00866 bufsize = strlen(neumCubeString);
00867 Vnm_print(2, "Vfetk_genCube: WARNING! USING EXPERIMENTAL NEUMANN BOUNDARY CONDITIONS!\n");
00868 VASSERT( bufsize <= VMAX_BUFSIZE );
00869 strncpy(buf, neumCubeString, VMAX_BUFSIZE);
00870 break;
00871 case VML_EXTERNAL:
00872 Vnm_print(2, "Vfetk_genCube: Got request for external mesh!\n");
00873 Vnm_print(2, "Vfetk_genCube: How did we get here?\n");
00874 return VRC_FAILURE;
00875 break;
00876 default:
00877 Vnm_print(2, "Vfetk_genCube: Unknown mesh type (%d)\n", meshType);
00878 return VRC_FAILURE;
00879 }
00880 VASSERT( VNULL != (sock=Vio_socketOpen(key,iodev,iofmt,iohost,iofile)) );
00881 Vio_bufTake(sock, buf, bufsize);
00882 AM_read(am, skey, sock);
00883 Vio_connectFree(sock);
00884 Vio_bufGive(sock);
00885 Vio_dtor(&sock);
00886
00887
00888 for (i=0; i<Gem_numVV(gm); i++) {
00889 vx = Gem_VV(gm, i);
00890 for (j=0; j<3; j++) {
00891 x = VV_coord(vx, j);
00892 x *= length[j];
00893 VV_setCoord(vx, j, x);
00894 }
00895 }
00896
00897
00898 for (i=0; i<Gem_numVV(gm); i++) {
00899 vx = Gem_VV(gm, i);
00900 for (j=0; j<3; j++) {
00901 x = VV_coord(vx, j);
00902 x += center[j];
00903 VV_setCoord(vx, j, x);
00904 }
00905 }
00906
00907
00908 return VRC_SUCCESS;
00909 }
00910
00911 VPUBLIC Vrc_Codes Vfetk_loadMesh(
00912 Vfetk *thee,
00913 double center[3],
00914 double length[3],
00915 Vfetk_MeshLoad meshType,
00916 Vio *sock) {
00917
00918 Vrc_Codes vrc;
00919 int skey = 0;
00920
00921
00922 switch (meshType) {
00923 case VML_EXTERNAL:
00924 if (sock == VNULL) {
00925 Vnm_print(2, "Vfetk_loadMesh: Got NULL socket!\n");
00926 return VRC_FAILURE;
00927 }
00928 AM_read(thee->am, skey, sock);
00929 Vio_connectFree(sock);
00930 Vio_bufGive(sock);
00931 Vio_dtor(&sock);
00932 break;
00933 case VML_DIRICUBE:
00934 vrc = Vfetk_genCube(thee, center, length, meshType);
00935 if (vrc == VRC_FAILURE) return VRC_FAILURE;
00936 break;
00937 case VML_NEUMCUBE:
00938 vrc = Vfetk_genCube(thee, center, length, meshType);
00939 if (vrc == VRC_FAILURE) return VRC_FAILURE;
00940 break;
00941 default:
00942 Vnm_print(2, "Vfetk_loadMesh: unrecognized mesh type (%d)!\n",
00943 meshType);
00944 return VRC_FAILURE;
00945 };
00946
00947
00948 Vnm_print(0, "Vfetk_ctor2: Constructing Vcsm...\n");
00949 thee->csm = VNULL;
00950 thee->csm = Vcsm_ctor(Vpbe_getValist(thee->pbe), thee->gm);
00951 VASSERT(thee->csm != VNULL);
00952 Vcsm_init(thee->csm);
00953
00954 return VRC_SUCCESS;
00955 }
00956
00957
00958 VPUBLIC void Bmat_printHB( Bmat *thee, char *fname ) {
00959
00960 Mat *Ablock;
00961 MATsym pqsym;
00962 int i, j, jj;
00963 int *IA, *JA;
00964 double *D, *L, *U;
00965 FILE *fp;
00966
00967 char mmtitle[72];
00968 char mmkey[] = {"8charkey"};
00969 int totc = 0, ptrc = 0, indc = 0, valc = 0;
00970 char mxtyp[] = {"RUA"};
00971 int nrow = 0, ncol = 0, numZ = 0;
00972 int numZdigits = 0, nrowdigits = 0;
00973 int nptrline = 8, nindline = 8, nvalline = 5;
00974 char ptrfmt[] = {"(8I10) "}, ptrfmtstr[] = {"%10d"};
00975 char indfmt[] = {"(8I10) "}, indfmtstr[] = {"%10d"};
00976 char valfmt[] = {"(5E16.8) "}, valfmtstr[] = {"%16.8E"};
00977
00978 VASSERT( thee->numB == 1 );
00979 Ablock = thee->AD[0][0];
00980
00981 VASSERT( Mat_format( Ablock ) == DRC_FORMAT );
00982
00983 pqsym = Mat_sym( Ablock );
00984
00985 if ( pqsym == IS_SYM ) {
00986 mxtyp[1] = 'S';
00987 } else if ( pqsym == ISNOT_SYM ) {
00988 mxtyp[1] = 'U';
00989 } else {
00990 VASSERT( 0 );
00991 }
00992
00993 nrow = Bmat_numRT( thee );
00994 ncol = Bmat_numCT( thee );
00995 numZ = Bmat_numZT( thee );
00996
00997 nrowdigits = (int) (log( nrow )/log( 10 )) + 1;
00998 numZdigits = (int) (log( numZ )/log( 10 )) + 1;
00999
01000 nptrline = (int) ( 80 / (numZdigits + 1) );
01001 nindline = (int) ( 80 / (nrowdigits + 1) );
01002
01003 sprintf(ptrfmt,"(%dI%d)",nptrline,numZdigits+1);
01004 sprintf(ptrfmtstr,"%%%dd",numZdigits+1);
01005 sprintf(indfmt,"(%dI%d)",nindline,nrowdigits+1);
01006 sprintf(indfmtstr,"%%%dd",nrowdigits+1);
01007
01008 ptrc = (int) ( ( (ncol + 1) - 1 ) / nptrline ) + 1;
01009 indc = (int) ( (numZ - 1) / nindline ) + 1;
01010 valc = (int) ( (numZ - 1) / nvalline ) + 1;
01011
01012 totc = ptrc + indc + valc;
01013
01014 sprintf( mmtitle, "Sparse '%s' Matrix - Harwell-Boeing Format - '%s'",
01015 thee->name, fname );
01016
01017
01018
01019 fp = fopen( fname, "w" );
01020 if (fp == VNULL) {
01021 Vnm_print(2,"Bmat_printHB: Ouch couldn't open file <%s>\n",fname);
01022 return;
01023 }
01024
01025
01026
01027 fprintf( fp, "%-72s%-8s\n", mmtitle, mmkey );
01028 fprintf( fp, "%14d%14d%14d%14d%14d\n", totc, ptrc, indc, valc, 0 );
01029 fprintf( fp, "%3s%11s%14d%14d%14d\n", mxtyp, " ", nrow, ncol, numZ );
01030 fprintf( fp, "%-16s%-16s%-20s%-20s\n", ptrfmt, indfmt, valfmt, "6E13.5" );
01031
01032 IA = Ablock->IA;
01033 JA = Ablock->JA;
01034 D = Ablock->diag;
01035 L = Ablock->offL;
01036 U = Ablock->offU;
01037
01038 if ( pqsym == IS_SYM ) {
01039
01040
01041
01042 for (i=0; i<(ncol+1); i++) {
01043 fprintf( fp, ptrfmtstr, Ablock->IA[i] + (i+1) );
01044 if ( ( (i+1) % nptrline ) == 0 ) {
01045 fprintf( fp, "\n" );
01046 }
01047 }
01048
01049 if ( ( (ncol+1) % nptrline ) != 0 ) {
01050 fprintf( fp, "\n" );
01051 }
01052
01053
01054
01055 j = 0;
01056 for (i=0; i<ncol; i++) {
01057 fprintf( fp, indfmtstr, i+1);
01058 if ( ( (j+1) % nindline ) == 0 ) {
01059 fprintf( fp, "\n" );
01060 }
01061 j++;
01062 for (jj=IA[i]; jj<IA[i+1]; jj++) {
01063 fprintf( fp, indfmtstr, JA[jj] + 1 );
01064 if ( ( (j+1) % nindline ) == 0 ) {
01065 fprintf( fp, "\n" );
01066 }
01067 j++;
01068 }
01069 }
01070
01071 if ( ( j % nindline ) != 0 ) {
01072 fprintf( fp, "\n" );
01073 }
01074
01075
01076
01077 j = 0;
01078 for (i=0; i<ncol; i++) {
01079 fprintf( fp, valfmtstr, D[i] );
01080 if ( ( (j+1) % nvalline ) == 0 ) {
01081 fprintf( fp, "\n" );
01082 }
01083 j++;
01084 for (jj=IA[i]; jj<IA[i+1]; jj++) {
01085 fprintf( fp, valfmtstr, L[jj] );
01086 if ( ( (j+1) % nvalline ) == 0 ) {
01087 fprintf( fp, "\n" );
01088 }
01089 j++;
01090 }
01091 }
01092
01093 if ( ( j % nvalline ) != 0 ) {
01094 fprintf( fp, "\n" );
01095 }
01096
01097 } else {
01098
01099 VASSERT( 0 );
01100 }
01101
01102
01103 fclose( fp );
01104 }
01105
01106 VPUBLIC PDE* Vfetk_PDE_ctor(Vfetk *fetk) {
01107
01108 PDE *thee = VNULL;
01109
01110 thee = Vmem_malloc(fetk->vmem, 1, sizeof(PDE));
01111 VASSERT(thee != VNULL);
01112 VASSERT(Vfetk_PDE_ctor2(thee, fetk));
01113
01114 return thee;
01115 }
01116
01117 VPUBLIC int Vfetk_PDE_ctor2(PDE *thee, Vfetk *fetk) {
01118
01119 int i;
01120
01121 if (thee == VNULL) {
01122 Vnm_print(2, "Vfetk_PDE_ctor2: Got NULL thee!\n");
01123 return 0;
01124 }
01125
01126
01127 var.fetk = fetk;
01128
01129
01130 thee->initAssemble = Vfetk_PDE_initAssemble;
01131 thee->initElement = Vfetk_PDE_initElement;
01132 thee->initFace = Vfetk_PDE_initFace;
01133 thee->initPoint = Vfetk_PDE_initPoint;
01134 thee->Fu = Vfetk_PDE_Fu;
01135 thee->Fu_v = Vfetk_PDE_Fu_v;
01136 thee->DFu_wv = Vfetk_PDE_DFu_wv;
01137 thee->delta = Vfetk_PDE_delta;
01138 thee->u_D = Vfetk_PDE_u_D;
01139 thee->u_T = Vfetk_PDE_u_T;
01140 thee->Ju = Vfetk_PDE_Ju;
01141 thee->vec = 1;
01142 thee->sym[0][0] = 1;
01143 thee->est[0] = 1.0;
01144 for (i=0; i<VMAX_BDTYPE; i++) thee->bmap[0][i] = i;
01145
01146
01147 thee->bisectEdge = Vfetk_PDE_bisectEdge;
01148 thee->mapBoundary = Vfetk_PDE_mapBoundary;
01149 thee->markSimplex = Vfetk_PDE_markSimplex;
01150 thee->oneChart = Vfetk_PDE_oneChart;
01151
01152
01153 thee->simplexBasisInit = Vfetk_PDE_simplexBasisInit;
01154 thee->simplexBasisForm = Vfetk_PDE_simplexBasisForm;
01155
01156 return 1;
01157 }
01158
01159 VPUBLIC void Vfetk_PDE_dtor(PDE **thee) {
01160
01161 if ((*thee) != VNULL) {
01162 Vfetk_PDE_dtor2(*thee);
01163
01164
01165
01166
01167
01168
01169 (*thee) = VNULL;
01170 }
01171
01172 }
01173
01174 VPUBLIC void Vfetk_PDE_dtor2(PDE *thee) {
01175 var.fetk = VNULL;
01176 }
01177
01178 VPRIVATE double smooth(int nverts, double dist[VAPBS_NVS], double coeff[VAPBS_NVS], int meth) {
01179
01180 int i;
01181 double weight;
01182 double num = 0.0;
01183 double den = 0.0;
01184
01185 for (i=0; i<nverts; i++) {
01186 if (dist[i] < VSMALL) return coeff[i];
01187 weight = 1.0/dist[i];
01188 if (meth == 0) {
01189 num += (weight * coeff[i]);
01190 den += weight;
01191 } else if (meth == 1) {
01192
01193
01194 if (coeff[i] < VSMALL) {
01195 num = 0.0;
01196 break;
01197 } else {
01198 num += weight; den += (weight/coeff[i]);
01199 }
01200 } else VASSERT(0);
01201 }
01202
01203 return (num/den);
01204
01205 }
01206
01207 VPRIVATE double diel() {
01208
01209 int i, j;
01210 double eps, epsp, epsw, dist[5], coeff[5], srad, swin, *vx;
01211 Vsurf_Meth srfm;
01212 Vacc *acc;
01213 PBEparm *pbeparm;
01214
01215 epsp = Vpbe_getSoluteDiel(var.fetk->pbe);
01216 epsw = Vpbe_getSolventDiel(var.fetk->pbe);
01217 VASSERT(var.fetk->pbeparm != VNULL);
01218 pbeparm = var.fetk->pbeparm;
01219 srfm = pbeparm->srfm;
01220 srad = pbeparm->srad;
01221 swin = pbeparm->swin;
01222 acc = var.fetk->pbe->acc;
01223
01224 eps = 0;
01225
01226 if (VABS(epsp - epsw) < VSMALL) return epsp;
01227 switch (srfm) {
01228 case VSM_MOL:
01229 eps = ((epsw-epsp)*Vacc_molAcc(acc, var.xq, srad) + epsp);
01230 break;
01231 case VSM_MOLSMOOTH:
01232 for (i=0; i<var.nverts; i++) {
01233 dist[i] = 0;
01234 vx = var.vx[i];
01235 for (j=0; j<3; j++) {
01236 dist[i] += VSQR(var.xq[j] - vx[j]);
01237 }
01238 dist[i] = VSQRT(dist[i]);
01239 coeff[i] = (epsw-epsp)*Vacc_molAcc(acc, var.xq, srad) + epsp;
01240 }
01241 eps = smooth(var.nverts, dist, coeff, 1);
01242 break;
01243 case VSM_SPLINE:
01244 eps = ((epsw-epsp)*Vacc_splineAcc(acc, var.xq, swin, 0.0) + epsp);
01245 break;
01246 default:
01247 Vnm_print(2, "Undefined surface method (%d)!\n", srfm);
01248 VASSERT(0);
01249 }
01250
01251 return eps;
01252 }
01253
01254 VPRIVATE double ionacc() {
01255
01256 int i, j;
01257 double dist[5], coeff[5], irad, swin, *vx, accval;
01258 Vsurf_Meth srfm;
01259 Vacc *acc = VNULL;
01260 PBEparm *pbeparm = VNULL;
01261
01262 VASSERT(var.fetk->pbeparm != VNULL);
01263 pbeparm = var.fetk->pbeparm;
01264 srfm = pbeparm->srfm;
01265 irad = Vpbe_getMaxIonRadius(var.fetk->pbe);
01266 swin = pbeparm->swin;
01267 acc = var.fetk->pbe->acc;
01268
01269 if (var.zks2 < VSMALL) return 0.0;
01270 switch (srfm) {
01271 case VSM_MOL:
01272 accval = Vacc_ivdwAcc(acc, var.xq, irad);
01273 break;
01274 case VSM_MOLSMOOTH:
01275 for (i=0; i<var.nverts; i++) {
01276 dist[i] = 0;
01277 vx = var.vx[i];
01278 for (j=0; j<3; j++) {
01279 dist[i] += VSQR(var.xq[j] - vx[j]);
01280 }
01281 dist[i] = VSQRT(dist[i]);
01282 coeff[i] = Vacc_ivdwAcc(acc, var.xq, irad);
01283 }
01284 accval = smooth(var.nverts, dist, coeff, 1);
01285 break;
01286 case VSM_SPLINE:
01287 accval = Vacc_splineAcc(acc, var.xq, swin, irad);
01288 break;
01289 default:
01290 Vnm_print(2, "Undefined surface method (%d)!\n", srfm);
01291 VASSERT(0);
01292 }
01293
01294 return accval;
01295 }
01296
01297 VPRIVATE double debye_U(Vpbe *pbe, int d, double x[]) {
01298
01299 double size, *position, charge, xkappa, eps_w, dist, T, pot, val;
01300 int iatom, i;
01301 Valist *alist;
01302 Vatom *atom;
01303
01304 eps_w = Vpbe_getSolventDiel(pbe);
01305 xkappa = (1.0e10)*Vpbe_getXkappa(pbe);
01306 T = Vpbe_getTemperature(pbe);
01307 alist = Vpbe_getValist(pbe);
01308 val = 0;
01309 pot = 0;
01310
01311 for (iatom=0; iatom<Valist_getNumberAtoms(alist); iatom++) {
01312 atom = Valist_getAtom(alist, iatom);
01313 position = Vatom_getPosition(atom);
01314 charge = Vunit_ec*Vatom_getCharge(atom);
01315 size = (1e-10)*Vatom_getRadius(atom);
01316 dist = 0;
01317 for (i=0; i<d; i++) {
01318 dist += VSQR(position[i] - x[i]);
01319 }
01320 dist = (1.0e-10)*VSQRT(dist);
01321 val = (charge)/(4*VPI*Vunit_eps0*eps_w*dist);
01322 if (xkappa != 0.0) {
01323 val = val*(exp(-xkappa*(dist-size))/(1+xkappa*size));
01324 }
01325 pot = pot + val;
01326 }
01327 pot = pot*Vunit_ec/(Vunit_kb*T);
01328
01329 return pot;
01330 }
01331
01332 VPRIVATE double debye_Udiff(Vpbe *pbe, int d, double x[]) {
01333
01334 double size, *position, charge, eps_p, dist, T, pot, val;
01335 double Ufull;
01336 int iatom, i;
01337 Valist *alist;
01338 Vatom *atom;
01339
01340 Ufull = debye_U(pbe, d, x);
01341
01342 eps_p = Vpbe_getSoluteDiel(pbe);
01343 T = Vpbe_getTemperature(pbe);
01344 alist = Vpbe_getValist(pbe);
01345 val = 0;
01346 pot = 0;
01347
01348 for (iatom=0; iatom<Valist_getNumberAtoms(alist); iatom++) {
01349 atom = Valist_getAtom(alist, iatom);
01350 position = Vatom_getPosition(atom);
01351 charge = Vunit_ec*Vatom_getCharge(atom);
01352 size = (1e-10)*Vatom_getRadius(atom);
01353 dist = 0;
01354 for (i=0; i<d; i++) {
01355 dist += VSQR(position[i] - x[i]);
01356 }
01357 dist = (1.0e-10)*VSQRT(dist);
01358 val = (charge)/(4*VPI*Vunit_eps0*eps_p*dist);
01359 pot = pot + val;
01360 }
01361 pot = pot*Vunit_ec/(Vunit_kb*T);
01362
01363 pot = Ufull - pot;
01364
01365 return pot;
01366 }
01367
01368 VPRIVATE void coulomb(Vpbe *pbe, int d, double pt[], double eps, double *U,
01369 double dU[], double *d2U) {
01370
01371 int iatom, i;
01372 double T, pot, fx, fy, fz, x, y, z, scale;
01373 double *position, charge, dist, dist2, val, vec[3], dUold[3], Uold;
01374 Valist *alist;
01375 Vatom *atom;
01376
01377
01378 T = Vpbe_getTemperature(pbe);
01379 alist = Vpbe_getValist(pbe);
01380 pot = 0; fx = 0; fy = 0; fz = 0;
01381 x = pt[0]; y = pt[1]; z = pt[2];
01382
01383
01384 if (!Vgreen_coulombD(var.green, 1, &x, &y, &z, &pot, &fx, &fy, &fz)) {
01385 Vnm_print(2, "Error calculating Green's function!\n");
01386 VASSERT(0);
01387 }
01388
01389
01390
01391 scale = Vunit_ec/(eps*Vunit_kb*T);
01392 *U = pot*scale;
01393 *d2U = 0.0;
01394 dU[0] = -fx*scale;
01395 dU[1] = -fy*scale;
01396 dU[2] = -fz*scale;
01397
01398 #if 0
01399
01400 val = 0.0;
01401 Uold = 0.0; dUold[0] = 0.0; dUold[1] = 0.0; dUold[2] = 0.0;
01402 for (iatom=0; iatom<Valist_getNumberAtoms(alist); iatom++) {
01403 atom = Valist_getAtom(alist, iatom);
01404 position = Vatom_getPosition(atom);
01405 charge = Vatom_getCharge(atom);
01406 dist2 = 0;
01407 for (i=0; i<d; i++) {
01408 vec[i] = (position[i] - pt[i]);
01409 dist2 += VSQR(vec[i]);
01410 }
01411 dist = VSQRT(dist2);
01412
01413
01414 Uold = Uold + charge/dist;
01415
01416
01417 for (i=0; i<d; i++) dUold[i] = dUold[i] + vec[i]*charge/(dist2*dist);
01418
01419 }
01420 Uold = Uold*VSQR(Vunit_ec)*(1.0e10)/(4*VPI*Vunit_eps0*eps*Vunit_kb*T);
01421 for (i=0; i<d; i++) {
01422 dUold[i] = dUold[i]*VSQR(Vunit_ec)*(1.0e10)/(4*VPI*Vunit_eps0*eps*Vunit_kb*T);
01423 }
01424
01425 printf("Unew - Uold = %g - %g = %g\n", *U, Uold, (*U - Uold));
01426 printf("||dUnew - dUold||^2 = %g\n", (VSQR(dU[0] - dUold[0])
01427 + VSQR(dU[1] - dUold[1]) + VSQR(dU[2] - dUold[2])));
01428 printf("dUnew[0] = %g, dUold[0] = %g\n", dU[0], dUold[0]);
01429 printf("dUnew[1] = %g, dUold[1] = %g\n", dU[1], dUold[1]);
01430 printf("dUnew[2] = %g, dUold[2] = %g\n", dU[2], dUold[2]);
01431
01432 #endif
01433
01434 }
01435
01436 VPUBLIC void Vfetk_PDE_initAssemble(PDE *thee, int ip[], double rp[]) {
01437
01438 #if 1
01439
01440
01441 if (var.initGreen) {
01442 Vgreen_dtor(&(var.green));
01443 var.initGreen = 0;
01444 }
01445 var.green = Vgreen_ctor(var.fetk->pbe->alist);
01446 var.initGreen = 1;
01447 #else
01448 if (!var.initGreen) {
01449 var.green = Vgreen_ctor(var.fetk->pbe->alist);
01450 var.initGreen = 1;
01451 }
01452 #endif
01453
01454 }
01455
01456 VPUBLIC void Vfetk_PDE_initElement(PDE *thee, int elementType, int chart,
01457 double tvx[][3], void *data) {
01458
01459 int i, j;
01460 double epsp, epsw;
01461
01462
01463
01464 VASSERT(data != NULL);
01465 var.simp = (SS *)data;
01466
01467
01468 var.sType = elementType;
01469
01470
01471 var.nverts = thee->dim+1;
01472 for (i=0; i<thee->dim+1; i++) var.verts[i] = SS_vertex(var.simp, i);
01473
01474
01475 for (i=0; i<thee->dim+1; i++) {
01476 for (j=0; j<thee->dim; j++) {
01477 var.vx[i][j] = tvx[i][j];
01478 }
01479 }
01480
01481
01482
01483
01484 var.jumpDiel = 0;
01485 }
01486
01487 VPUBLIC void Vfetk_PDE_initFace(PDE *thee, int faceType, int chart,
01488 double tnvec[]) {
01489
01490 int i;
01491
01492
01493 for (i=0; i<thee->dim; i++) var.nvec[i] = tnvec[i];
01494
01495
01496 var.fType = faceType;
01497 }
01498
01499 VPUBLIC void Vfetk_PDE_initPoint(PDE *thee, int pointType, int chart,
01500 double txq[], double tU[], double tdU[][3]) {
01501
01502 int i, j, ichop;
01503 double u2, coef2, eps_p;
01504 Vhal_PBEType pdetype;
01505 Vpbe *pbe = VNULL;
01506
01507 eps_p = Vpbe_getSoluteDiel(var.fetk->pbe);
01508 pdetype = var.fetk->type;
01509 pbe = var.fetk->pbe;
01510
01511
01512
01513 if ((pdetype == PBE_LRPBE) || (pdetype == PBE_NRPBE)) {
01514 coulomb(pbe, thee->dim, txq, eps_p, &(var.W), var.dW, &(var.d2W));
01515 }
01516 for (i=0; i<thee->vec; i++) {
01517 var.U[i] = tU[i];
01518 for (j=0; j<thee->dim; j++) {
01519 var.xq[j] = txq[j];
01520 var.dU[i][j] = tdU[i][j];
01521 }
01522 }
01523
01524
01525 if (pointType == 0) {
01526
01527
01528 var.diel = diel();
01529 var.ionacc = ionacc();
01530 var.A = var.diel;
01531 var.F = (var.diel - eps_p);
01532
01533 switch (pdetype) {
01534
01535 case PBE_LPBE:
01536 var.DB = var.ionacc*var.zkappa2*var.ionstr;
01537 var.B = var.DB*var.U[0];
01538 break;
01539
01540 case PBE_NPBE:
01541
01542 var.B = 0;
01543 var.DB = 0;
01544 if ((var.ionacc > VSMALL) && (var.zks2 > VSMALL)) {
01545 for (i=0; i<var.nion; i++) {
01546 u2 = -1.0 * var.U[0] * var.ionQ[i];
01547
01548
01549 coef2 = -1.0 * var.ionacc * var.zks2 * var.ionConc[i];
01550 var.B += (coef2 * Vcap_exp(u2, &ichop));
01551
01552 coef2 = -1.0 * var.ionQ[i] * coef2;
01553 var.DB += (coef2 * Vcap_exp(u2, &ichop));
01554 }
01555 }
01556 break;
01557
01558 case PBE_LRPBE:
01559 var.DB = var.ionacc*var.zkappa2*var.ionstr;
01560 var.B = var.DB*(var.U[0]+var.W);
01561 break;
01562
01563 case PBE_NRPBE:
01564
01565 var.B = 0;
01566 var.DB = 0;
01567 if ((var.ionacc > VSMALL) && (var.zks2 > VSMALL)) {
01568 for (i=0; i<var.nion; i++) {
01569 u2 = -1.0 * (var.U[0] + var.W) * var.ionQ[i];
01570
01571
01572 coef2 = -1.0 * var.ionacc * var.zks2 * var.ionConc[i];
01573 var.B += (coef2 * Vcap_exp(u2, &ichop));
01574
01575
01576 coef2 = -1.0 * var.ionQ[i] * coef2;
01577 var.DB += (coef2 * Vcap_exp(u2, &ichop));
01578 }
01579 }
01580 break;
01581
01582 case PBE_SMPBE:
01583
01584 var.B = 0;
01585 var.DB = 0;
01586 if ((var.ionacc > VSMALL) && (var.zks2 > VSMALL)) {
01587 for (i=0; i<var.nion; i++) {
01588 u2 = -1.0 * var.U[0] * var.ionQ[i];
01589
01590
01591 coef2 = -1.0 * var.ionacc * var.zks2 * var.ionConc[i];
01592 var.B += (coef2 * Vcap_exp(u2, &ichop));
01593
01594 coef2 = -1.0 * var.ionQ[i] * coef2;
01595 var.DB += (coef2 * Vcap_exp(u2, &ichop));
01596 }
01597 }
01598 break;
01599 default:
01600 Vnm_print(2, "Vfetk_PDE_initPoint: Unknown PBE type (%d)!\n",
01601 pdetype);
01602 VASSERT(0);
01603 break;
01604 }
01605
01606
01607
01608 } else {
01609 #ifdef DONEUMANN
01610 ;
01611 #else
01612 Vnm_print(2, "Vfetk: Whoa! I just got a boundary point to evaluate (%d)!\n", pointType);
01613 Vnm_print(2, "Vfetk: Did you do that on purpose?\n");
01614 #endif
01615 }
01616
01617 #if 0
01618 Vfetk_dumpLocalVar();
01619 #endif
01620
01621 }
01622
01623 VPUBLIC void Vfetk_PDE_Fu(PDE *thee, int key, double F[]) {
01624
01625
01626
01627 F[0] = 0.;
01628
01629 }
01630
01631 VPUBLIC double Vfetk_PDE_Fu_v(
01632 PDE *thee,
01633 int key,
01634 double V[],
01635 double dV[][VAPBS_DIM]
01636 ) {
01637
01638 Vhal_PBEType type;
01639 int i;
01640 double value = 0.;
01641
01642 type = var.fetk->type;
01643
01644
01645 if (key == 0) {
01646
01647 for (i=0; i<thee->dim; i++) value += ( var.A * var.dU[0][i] * dV[0][i] );
01648 value += var.B * V[0];
01649
01650 if ((type == PBE_LRPBE) || (type == PBE_NRPBE)) {
01651 for (i=0; i<thee->dim; i++) {
01652 if (var.F > VSMALL) value += (var.F * var.dW[i] * dV[0][i]);
01653 }
01654 }
01655
01656
01657 } else {
01658 #ifdef DONEUMANN
01659 value = 0.0;
01660 #else
01661 Vnm_print(2, "Vfetk: Whoa! I was just asked to evaluate a boundary weak form for point type %d!\n", key);
01662 #endif
01663 }
01664
01665 var.Fu_v = value;
01666 return value;
01667 }
01668
01669 VPUBLIC double Vfetk_PDE_DFu_wv(
01670 PDE *thee,
01671 int key,
01672 double W[],
01673 double dW[][VAPBS_DIM],
01674 double V[],
01675 double dV[][3]
01676 ) {
01677
01678 Vhal_PBEType type;
01679 int i;
01680 double value = 0.;
01681
01682 type = var.fetk->type;
01683
01684
01685 if (key == 0) {
01686 value += var.DB * W[0] * V[0];
01687 for (i=0; i<thee->dim; i++) value += ( var.A * dW[0][i] * dV[0][i] );
01688
01689
01690 } else {
01691 #ifdef DONEUMANN
01692 value = 0.0;
01693 #else
01694 Vnm_print(2, "Vfetk: Whoa! I was just asked to evaluate a boundary weak form for point type %d!\n", key);
01695 #endif
01696 }
01697
01698 var.DFu_wv = value;
01699 return value;
01700 }
01701
01704 #define VRINGMAX 1000
01705
01707 #define VATOMMAX 1000000
01708 VPUBLIC void Vfetk_PDE_delta(PDE *thee, int type, int chart, double txq[],
01709 void *user, double F[]) {
01710
01711 int iatom, jatom, natoms, atomIndex, atomList[VATOMMAX], nAtomList;
01712 int gotAtom, numSring, isimp, ivert, sid;
01713 double *position, charge, phi[VAPBS_NVS], phix[VAPBS_NVS][3], value;
01714 Vatom *atom;
01715 Vhal_PBEType pdetype;
01716 SS *sring[VRINGMAX];
01717 VV *vertex = (VV *)user;
01718
01719 pdetype = var.fetk->type;
01720
01721 F[0] = 0.0;
01722
01723 if ((pdetype == PBE_LPBE) || (pdetype == PBE_NPBE) || (pdetype == PBE_SMPBE) ) {
01724 VASSERT( vertex != VNULL);
01725 numSring = 0;
01726 sring[numSring] = VV_firstSS(vertex);
01727 while (sring[numSring] != VNULL) {
01728 numSring++;
01729 sring[numSring] = SS_link(sring[numSring-1], vertex);
01730 }
01731 VASSERT( numSring > 0 );
01732 VASSERT( numSring <= VRINGMAX );
01733
01734
01735 F[0] = 0.;
01736 charge = 0.;
01737 nAtomList = 0;
01738 for (isimp=0; isimp<numSring; isimp++) {
01739 sid = SS_id(sring[isimp]);
01740 natoms = Vcsm_getNumberAtoms(Vfetk_getVcsm(var.fetk), sid);
01741 for (iatom=0; iatom<natoms; iatom++) {
01742
01743 atomIndex = Vcsm_getAtomIndex(Vfetk_getVcsm(var.fetk),
01744 iatom, sid);
01745 gotAtom = 0;
01746 for (jatom=0; jatom<nAtomList; jatom++) {
01747 if (atomList[jatom] == atomIndex) {
01748 gotAtom = 1;
01749 break;
01750 }
01751 }
01752 if (!gotAtom) {
01753 VASSERT(nAtomList < VATOMMAX);
01754 atomList[nAtomList] = atomIndex;
01755 nAtomList++;
01756
01757 atom = Vcsm_getAtom(Vfetk_getVcsm(var.fetk), iatom, sid);
01758 charge = Vatom_getCharge(atom);
01759 position = Vatom_getPosition(atom);
01760
01761
01762
01763
01764
01765 if (!Gem_pointInSimplexVal(Vfetk_getGem(var.fetk),
01766 sring[isimp], position, phi, phix)) {
01767 if (!Gem_pointInSimplex(Vfetk_getGem(var.fetk),
01768 sring[isimp], position)) {
01769 Vnm_print(2, "delta: Both Gem_pointInSimplexVal \
01770 and Gem_pointInSimplex detected misplaced point charge!\n");
01771 Vnm_print(2, "delta: I think you have problems: \
01772 phi = {");
01773 for (ivert=0; ivert<Gem_dimVV(Vfetk_getGem(var.fetk)); ivert++) Vnm_print(2, "%e ", phi[ivert]);
01774 Vnm_print(2, "}\n");
01775 }
01776 }
01777 value = 0;
01778 for (ivert=0; ivert<Gem_dimVV(Vfetk_getGem(var.fetk)); ivert++) {
01779 if (VV_id(SS_vertex(sring[isimp], ivert)) == VV_id(vertex)) value += phi[ivert];
01780 }
01781
01782 F[0] += (value * Vpbe_getZmagic(var.fetk->pbe) * charge);
01783 }
01784 }
01785 }
01786
01787 } else if ((pdetype == PBE_LRPBE) || (pdetype == PBE_NRPBE)) {
01788 F[0] = 0.0;
01789 } else { VASSERT(0); }
01790
01791 var.delta = F[0];
01792
01793 }
01794
01795 VPUBLIC void Vfetk_PDE_u_D(PDE *thee, int type, int chart, double txq[],
01796 double F[]) {
01797
01798 if ((var.fetk->type == PBE_LPBE) || (var.fetk->type == PBE_NPBE) || (var.fetk->type == PBE_SMPBE) ) {
01799 F[0] = debye_U(var.fetk->pbe, thee->dim, txq);
01800 } else if ((var.fetk->type == PBE_LRPBE) || (var.fetk->type == PBE_NRPBE)) {
01801 F[0] = debye_Udiff(var.fetk->pbe, thee->dim, txq);
01802 } else VASSERT(0);
01803
01804 var.u_D = F[0];
01805
01806 }
01807
01808 VPUBLIC void Vfetk_PDE_u_T(PDE *thee, int type, int chart, double txq[],
01809 double F[]) {
01810
01811 F[0] = 0.0;
01812 var.u_T = F[0];
01813
01814 }
01815
01816
01817 VPUBLIC void Vfetk_PDE_bisectEdge(int dim, int dimII, int edgeType,
01818 int chart[], double vx[][3]) {
01819
01820 int i;
01821
01822 for (i=0; i<dimII; i++) vx[2][i] = .5 * (vx[0][i] + vx[1][i]);
01823 chart[2] = chart[0];
01824
01825 }
01826
01827 VPUBLIC void Vfetk_PDE_mapBoundary(int dim, int dimII, int vertexType,
01828 int chart, double vx[3]) {
01829
01830 }
01831
01832 VPUBLIC int Vfetk_PDE_markSimplex(int dim, int dimII, int simplexType,
01833 int faceType[VAPBS_NVS], int vertexType[VAPBS_NVS], int chart[], double vx[][3],
01834 void *simplex) {
01835
01836 double targetRes, edgeLength, srad, swin, myAcc, refAcc;
01837 int i, natoms;
01838 Vsurf_Meth srfm;
01839 Vhal_PBEType type;
01840 FEMparm *feparm = VNULL;
01841 PBEparm *pbeparm = VNULL;
01842 Vpbe *pbe = VNULL;
01843 Vacc *acc = VNULL;
01844 Vcsm *csm = VNULL;
01845 SS *simp = VNULL;
01846
01847 VASSERT(var.fetk->feparm != VNULL);
01848 feparm = var.fetk->feparm;
01849 VASSERT(var.fetk->pbeparm != VNULL);;
01850 pbeparm = var.fetk->pbeparm;
01851 pbe = var.fetk->pbe;
01852 csm = Vfetk_getVcsm(var.fetk);
01853 acc = pbe->acc;
01854 targetRes = feparm->targetRes;
01855 srfm = pbeparm->srfm;
01856 srad = pbeparm->srad;
01857 swin = pbeparm->swin;
01858 simp = (SS *)simplex;
01859 type = var.fetk->type;
01860
01861
01862
01863 Gem_longestEdge(var.fetk->gm, simp, -1, &edgeLength);
01864 if (edgeLength < targetRes) return 0;
01865
01866
01867 if ((type == PBE_LPBE) || (type == PBE_NPBE) || (type == PBE_SMPBE) ) {
01868 natoms = Vcsm_getNumberAtoms(csm, SS_id(simp));
01869 if (natoms > 0) {
01870 return 1;
01871 }
01872 }
01873
01874
01875
01876
01877 switch(srfm) {
01878 case VSM_MOL:
01879 refAcc = Vacc_molAcc(acc, vx[0], srad);
01880 for (i=1; i<(dim+1); i++) {
01881 myAcc = Vacc_molAcc(acc, vx[i], srad);
01882 if (myAcc != refAcc) {
01883 return 1;
01884 }
01885 }
01886 break;
01887 case VSM_MOLSMOOTH:
01888 refAcc = Vacc_molAcc(acc, vx[0], srad);
01889 for (i=1; i<(dim+1); i++) {
01890 myAcc = Vacc_molAcc(acc, vx[i], srad);
01891 if (myAcc != refAcc) {
01892 return 1;
01893 }
01894 }
01895 break;
01896 case VSM_SPLINE:
01897 refAcc = Vacc_splineAcc(acc, vx[0], swin, 0.0);
01898 for (i=1; i<(dim+1); i++) {
01899 myAcc = Vacc_splineAcc(acc, vx[i], swin, 0.0);
01900 if (myAcc != refAcc) {
01901 return 1;
01902 }
01903 }
01904 break;
01905 default:
01906 VASSERT(0);
01907 break;
01908 }
01909
01910 return 0;
01911 }
01912
01913 VPUBLIC void Vfetk_PDE_oneChart(int dim, int dimII, int objType, int chart[],
01914 double vx[][3], int dimV) {
01915
01916 }
01917
01918 VPUBLIC double Vfetk_PDE_Ju(PDE *thee, int key) {
01919
01920 int i, ichop;
01921 double dielE, qmE, coef2, u2;
01922 double value = 0.;
01923 Vhal_PBEType type;
01924
01925 type = var.fetk->type;
01926
01927
01928 if (key == 0) {
01929 dielE = 0;
01930 for (i=0; i<3; i++) dielE += VSQR(var.dU[0][i]);
01931 dielE = dielE*var.diel;
01932
01933 switch (type) {
01934 case PBE_LPBE:
01935 if ((var.ionacc > VSMALL) && (var.zkappa2 > VSMALL)) {
01936 qmE = var.ionacc*var.zkappa2*VSQR(var.U[0]);
01937 } else qmE = 0;
01938 break;
01939 case PBE_NPBE:
01940 if ((var.ionacc > VSMALL) && (var.zks2 > VSMALL)) {
01941 qmE = 0.;
01942 for (i=0; i<var.nion; i++) {
01943 coef2 = var.ionacc * var.zks2 * var.ionConc[i] * var.ionQ[i];
01944 u2 = -1.0 * (var.U[0]) * var.ionQ[i];
01945 qmE += (coef2 * (Vcap_exp(u2, &ichop) - 1.0));
01946 }
01947 } else qmE = 0;
01948 break;
01949 case PBE_LRPBE:
01950 if ((var.ionacc > VSMALL) && (var.zkappa2 > VSMALL)) {
01951 qmE = var.ionacc*var.zkappa2*VSQR((var.U[0] + var.W));
01952 } else qmE = 0;
01953 break;
01954 case PBE_NRPBE:
01955 if ((var.ionacc > VSMALL) && (var.zks2 > VSMALL)) {
01956 qmE = 0.;
01957 for (i=0; i<var.nion; i++) {
01958 coef2 = var.ionacc * var.zks2 * var.ionConc[i] * var.ionQ[i];
01959 u2 = -1.0 * (var.U[0] + var.W) * var.ionQ[i];
01960 qmE += (coef2 * (Vcap_exp(u2, &ichop) - 1.0));
01961 }
01962 } else qmE = 0;
01963 break;
01964 case PBE_SMPBE:
01965 if ((var.ionacc > VSMALL) && (var.zks2 > VSMALL)) {
01966 qmE = 0.;
01967 for (i=0; i<var.nion; i++) {
01968 coef2 = var.ionacc * var.zks2 * var.ionConc[i] * var.ionQ[i];
01969 u2 = -1.0 * (var.U[0]) * var.ionQ[i];
01970 qmE += (coef2 * (Vcap_exp(u2, &ichop) - 1.0));
01971 }
01972 } else qmE = 0;
01973 break;
01974 default:
01975 Vnm_print(2, "Vfetk_PDE_Ju: Invalid PBE type (%d)!\n", type);
01976 VASSERT(0);
01977 break;
01978 }
01979
01980 value = 0.5*(dielE + qmE)/Vpbe_getZmagic(var.fetk->pbe);
01981
01982
01983 } else if (key == 1) {
01984 value = 0.0;
01985
01986
01987 } else VASSERT(0);
01988
01989 return value;
01990
01991 }
01992
01993 VPUBLIC void Vfetk_externalUpdateFunction(SS **simps, int num) {
01994
01995 Vcsm *csm = VNULL;
01996 int rc;
01997
01998 VASSERT(var.fetk != VNULL);
01999 csm = Vfetk_getVcsm(var.fetk);
02000 VASSERT(csm != VNULL);
02001
02002 rc = Vcsm_update(csm, simps, num);
02003
02004 if (!rc) {
02005 Vnm_print(2, "Error while updating charge-simplex map!\n");
02006 VASSERT(0);
02007 }
02008 }
02009
02010 VPRIVATE void polyEval(int numP, double p[], double c[][VMAXP], double xv[]) {
02011 int i;
02012 double x, y, z;
02013
02014 x = xv[0];
02015 y = xv[1];
02016 z = xv[2];
02017 for (i=0; i<numP; i++) {
02018 p[i] = c[i][0]
02019 + c[i][1] * x
02020 + c[i][2] * y
02021 + c[i][3] * z
02022 + c[i][4] * x*x
02023 + c[i][5] * y*y
02024 + c[i][6] * z*z
02025 + c[i][7] * x*y
02026 + c[i][8] * x*z
02027 + c[i][9] * y*z
02028 + c[i][10] * x*x*x
02029 + c[i][11] * y*y*y
02030 + c[i][12] * z*z*z
02031 + c[i][13] * x*x*y
02032 + c[i][14] * x*x*z
02033 + c[i][15] * x*y*y
02034 + c[i][16] * y*y*z
02035 + c[i][17] * x*z*z
02036 + c[i][18] * y*z*z;
02037 }
02038 }
02039
02040 VPRIVATE void setCoef(int numP, double c[][VMAXP], double cx[][VMAXP],
02041 double cy[][VMAXP], double cz[][VMAXP], int ic[][VMAXP], int icx[][VMAXP],
02042 int icy[][VMAXP], int icz[][VMAXP]) {
02043
02044 int i, j;
02045 for (i=0; i<numP; i++) {
02046 for (j=0; j<VMAXP; j++) {
02047 c[i][j] = 0.5 * (double)ic[i][j];
02048 cx[i][j] = 0.5 * (double)icx[i][j];
02049 cy[i][j] = 0.5 * (double)icy[i][j];
02050 cz[i][j] = 0.5 * (double)icz[i][j];
02051 }
02052 }
02053 }
02054
02055 VPUBLIC int Vfetk_PDE_simplexBasisInit(int key, int dim, int comp, int *ndof,
02056 int dof[]) {
02057
02058 int qorder, bump, dimIS[VAPBS_NVS];
02059
02060
02061 qorder = P_DEG;
02062
02063
02064 if ((key == 0) || (key == 1)) {
02065 bump = 0;
02066 } else if ((key == 2) || (key == 3)) {
02067 bump = 1;
02068 } else { VASSERT(0); }
02069
02070
02071 if (dim==2) {
02072
02073 dimIS[0] = 3;
02074 dimIS[1] = 3;
02075 dimIS[2] = 0;
02076 dimIS[3] = 1;
02077 if (bump==0) {
02078 if (P_DEG==1) {
02079 init_2DP1(dimIS, ndof, dof, c, cx, cy, cz);
02080 } else if (P_DEG==2) {
02081 init_2DP1(dimIS, ndof, dof, c, cx, cy, cz);
02082 } else if (P_DEG==3) {
02083 init_2DP1(dimIS, ndof, dof, c, cx, cy, cz);
02084 } else Vnm_print(2, "..bad order..");
02085 } else if (bump==1) {
02086 if (P_DEG==1) {
02087 init_2DP1(dimIS, ndof, dof, c, cx, cy, cz);
02088 } else Vnm_print(2, "..bad order..");
02089 } else Vnm_print(2, "..bad bump..");
02090 } else if (dim==3) {
02091
02092 dimIS[0] = 4;
02093 dimIS[1] = 6;
02094 dimIS[2] = 4;
02095 dimIS[3] = 1;
02096 if (bump==0) {
02097 if (P_DEG==1) {
02098 init_3DP1(dimIS, ndof, dof, c, cx, cy, cz);
02099 } else if (P_DEG==2) {
02100 init_3DP1(dimIS, ndof, dof, c, cx, cy, cz);
02101 } else if (P_DEG==3) {
02102 init_3DP1(dimIS, ndof, dof, c, cx, cy, cz);
02103 } else Vnm_print(2, "..bad order..");
02104 } else if (bump==1) {
02105 if (P_DEG==1) {
02106 init_3DP1(dimIS, ndof, dof, c, cx, cy, cz);
02107 } else Vnm_print(2, "..bad order..");
02108 } else Vnm_print(2, "..bad bump..");
02109 } else Vnm_print(2, "..bad dimension..");
02110
02111
02112 numP = *ndof;
02113
02114
02115 return qorder;
02116 }
02117
02118 VPUBLIC void Vfetk_PDE_simplexBasisForm(int key, int dim, int comp, int pdkey,
02119 double xq[], double basis[]) {
02120
02121 if (pdkey == 0) {
02122 polyEval(numP, basis, c, xq);
02123 } else if (pdkey == 1) {
02124 polyEval(numP, basis, cx, xq);
02125 } else if (pdkey == 2) {
02126 polyEval(numP, basis, cy, xq);
02127 } else if (pdkey == 3) {
02128 polyEval(numP, basis, cz, xq);
02129 } else { VASSERT(0); }
02130 }
02131
02132 VPRIVATE void init_2DP1(int dimIS[], int *ndof, int dof[], double c[][VMAXP],
02133 double cx[][VMAXP], double cy[][VMAXP], double cz[][VMAXP]) {
02134
02135 int i;
02136
02137
02138 dof[0] = 1;
02139 dof[1] = 0;
02140 dof[2] = 0;
02141 dof[3] = 0;
02142 *ndof = 0;
02143 for (i=0; i<VAPBS_NVS; i++) *ndof += dimIS[i] * dof[i];
02144 VASSERT( *ndof == dim_2DP1 );
02145 VASSERT( *ndof <= VMAXP );
02146
02147
02148 setCoef( *ndof, c, cx, cy, cz, lgr_2DP1, lgr_2DP1x, lgr_2DP1y, lgr_2DP1z );
02149 }
02150
02151 VPRIVATE void init_3DP1(int dimIS[], int *ndof, int dof[], double c[][VMAXP],
02152 double cx[][VMAXP], double cy[][VMAXP], double cz[][VMAXP]) {
02153
02154 int i;
02155
02156
02157 dof[0] = 1;
02158 dof[1] = 0;
02159 dof[2] = 0;
02160 dof[3] = 0;
02161 *ndof = 0;
02162 for (i=0; i<VAPBS_NVS; i++) *ndof += dimIS[i] * dof[i];
02163 VASSERT( *ndof == dim_3DP1 );
02164 VASSERT( *ndof <= VMAXP );
02165
02166
02167 setCoef( *ndof, c, cx, cy, cz, lgr_3DP1, lgr_3DP1x, lgr_3DP1y, lgr_3DP1z );
02168 }
02169
02170 VPUBLIC void Vfetk_dumpLocalVar() {
02171
02172 int i;
02173
02174 Vnm_print(1, "DEBUG: nvec = (%g, %g, %g)\n", var.nvec[0], var.nvec[1],
02175 var.nvec[2]);
02176 Vnm_print(1, "DEBUG: nverts = %d\n", var.nverts);
02177 for (i=0; i<var.nverts; i++) {
02178 Vnm_print(1, "DEBUG: verts[%d] ID = %d\n", i, VV_id(var.verts[i]));
02179 Vnm_print(1, "DEBUG: vx[%d] = (%g, %g, %g)\n", i, var.vx[i][0],
02180 var.vx[i][1], var.vx[i][2]);
02181 }
02182 Vnm_print(1, "DEBUG: simp ID = %d\n", SS_id(var.simp));
02183 Vnm_print(1, "DEBUG: sType = %d\n", var.sType);
02184 Vnm_print(1, "DEBUG: fType = %d\n", var.fType);
02185 Vnm_print(1, "DEBUG: xq = (%g, %g, %g)\n", var.xq[0], var.xq[1], var.xq[2]);
02186 Vnm_print(1, "DEBUG: U[0] = %g\n", var.U[0]);
02187 Vnm_print(1, "DEBUG: dU[0] = (%g, %g, %g)\n", var.dU[0][0], var.dU[0][1],
02188 var.dU[0][2]);
02189 Vnm_print(1, "DEBUG: W = %g\n", var.W);
02190 Vnm_print(1, "DEBUG: d2W = %g\n", var.d2W);
02191 Vnm_print(1, "DEBUG: dW = (%g, %g, %g)\n", var.dW[0], var.dW[1], var.dW[2]);
02192 Vnm_print(1, "DEBUG: diel = %g\n", var.diel);
02193 Vnm_print(1, "DEBUG: ionacc = %g\n", var.ionacc);
02194 Vnm_print(1, "DEBUG: A = %g\n", var.A);
02195 Vnm_print(1, "DEBUG: F = %g\n", var.F);
02196 Vnm_print(1, "DEBUG: B = %g\n", var.B);
02197 Vnm_print(1, "DEBUG: DB = %g\n", var.DB);
02198 Vnm_print(1, "DEBUG: nion = %d\n", var.nion);
02199 for (i=0; i<var.nion; i++) {
02200 Vnm_print(1, "DEBUG: ionConc[%d] = %g\n", i, var.ionConc[i]);
02201 Vnm_print(1, "DEBUG: ionQ[%d] = %g\n", i, var.ionQ[i]);
02202 Vnm_print(1, "DEBUG: ionRadii[%d] = %g\n", i, var.ionRadii[i]);
02203 }
02204 Vnm_print(1, "DEBUG: zkappa2 = %g\n", var.zkappa2);
02205 Vnm_print(1, "DEBUG: zks2 = %g\n", var.zks2);
02206 Vnm_print(1, "DEBUG: Fu_v = %g\n", var.Fu_v);
02207 Vnm_print(1, "DEBUG: DFu_wv = %g\n", var.DFu_wv);
02208 Vnm_print(1, "DEBUG: delta = %g\n", var.delta);
02209 Vnm_print(1, "DEBUG: u_D = %g\n", var.u_D);
02210 Vnm_print(1, "DEBUG: u_T = %g\n", var.u_T);
02211
02212 };
02213
02214 VPUBLIC int Vfetk_fillArray(Vfetk *thee, Bvec *vec, Vdata_Type type) {
02215
02216 int i, j, ichop;
02217 double coord[3], chi, q, conc, val;
02218 VV *vert;
02219 Bvec *u, *u_d;
02220 AM *am;
02221 Gem *gm;
02222 PBEparm *pbeparm;
02223 Vacc *acc;
02224 Vpbe *pbe;
02225
02226 gm = thee->gm;
02227 am = thee->am;
02228 pbe = thee->pbe;
02229 pbeparm = thee->pbeparm;
02230 acc = pbe->acc;
02231
02232
02233 if (Bvec_numRB(vec, 0) != Gem_numVV(gm)) {
02234 Vnm_print(2, "Vfetk_fillArray: insufficient space in Bvec!\n");
02235 Vnm_print(2, "Vfetk_fillArray: Have %d, need %d!\n", Bvec_numRB(vec, 0),
02236 Gem_numVV(gm));
02237 return 0;
02238 }
02239
02240 switch (type) {
02241
02242 case VDT_CHARGE:
02243 Vnm_print(2, "Vfetk_fillArray: can't write out charge distribution!\n");
02244 return 0;
02245 break;
02246
02247 case VDT_POT:
02248 u = am->u;
02249 u_d = am->ud;
02250
02251 Bvec_copy(vec, u);
02252
02253 Bvec_axpy(vec, u_d, 1.0);
02254 break;
02255
02256 case VDT_SMOL:
02257 for (i=0; i<Gem_numVV(gm); i++) {
02258 vert = Gem_VV(gm, i);
02259 for (j=0; j<3; j++) coord[j] = VV_coord(vert, j);
02260 chi = Vacc_molAcc(acc, coord, pbe->solventRadius);
02261 Bvec_set(vec, i, chi);
02262 }
02263 break;
02264
02265 case VDT_SSPL:
02266 for (i=0; i<Gem_numVV(gm); i++) {
02267 vert = Gem_VV(gm, i);
02268 for (j=0; j<3; j++) coord[j] = VV_coord(vert, j);
02269 chi = Vacc_splineAcc(acc, coord, pbeparm->swin, 0.0);
02270 Bvec_set(vec, i, chi);
02271 }
02272 break;
02273
02274 case VDT_VDW:
02275 for (i=0; i<Gem_numVV(gm); i++) {
02276 vert = Gem_VV(gm, i);
02277 for (j=0; j<3; j++) coord[j] = VV_coord(vert, j);
02278 chi = Vacc_vdwAcc(acc, coord);
02279 Bvec_set(vec, i, chi);
02280 }
02281 break;
02282
02283 case VDT_IVDW:
02284 for (i=0; i<Gem_numVV(gm); i++) {
02285 vert = Gem_VV(gm, i);
02286 for (j=0; j<3; j++) coord[j] = VV_coord(vert, j);
02287 chi = Vacc_ivdwAcc(acc, coord, pbe->maxIonRadius);
02288 Bvec_set(vec, i, chi);
02289 }
02290 break;
02291
02292 case VDT_LAP:
02293 Vnm_print(2, "Vfetk_fillArray: can't write out Laplacian!\n");
02294 return 0;
02295 break;
02296
02297 case VDT_EDENS:
02298 Vnm_print(2, "Vfetk_fillArray: can't write out energy density!\n");
02299 return 0;
02300 break;
02301
02302 case VDT_NDENS:
02303 u = am->u;
02304 u_d = am->ud;
02305
02306 Bvec_copy(vec, u);
02307
02308 Bvec_axpy(vec, u_d, 1.0);
02309
02310 ichop = 0;
02311 for (i=0; i<Gem_numVV(gm); i++) {
02312 val = 0;
02313 for (j=0; j<pbe->numIon; j++) {
02314 q = pbe->ionQ[j];
02315 conc = pbe->ionConc[j];
02316 if (thee->type == PBE_NPBE || thee->type == PBE_SMPBE ) {
02317 val += (conc*Vcap_exp(-q*Bvec_val(vec, i), &ichop));
02318 } else if (thee->type == PBE_LPBE) {
02319 val += (conc * ( 1 - q*Bvec_val(vec, i)));
02320 }
02321 }
02322 Bvec_set(vec, i, val);
02323 }
02324 break;
02325
02326 case VDT_QDENS:
02327 u = am->u;
02328 u_d = am->ud;
02329
02330 Bvec_copy(vec, u);
02331
02332 Bvec_axpy(vec, u_d, 1.0);
02333
02334 ichop = 0;
02335 for (i=0; i<Gem_numVV(gm); i++) {
02336 val = 0;
02337 for (j=0; j<pbe->numIon; j++) {
02338 q = pbe->ionQ[j];
02339 conc = pbe->ionConc[j];
02340 if (thee->type == PBE_NPBE || thee->type == PBE_SMPBE ) {
02341 val += (q*conc*Vcap_exp(-q*Bvec_val(vec, i), &ichop));
02342 } else if (thee->type == PBE_LPBE) {
02343 val += (q*conc*(1 - q*Bvec_val(vec, i)));
02344 }
02345 }
02346 Bvec_set(vec, i, val);
02347 }
02348 break;
02349
02350 case VDT_DIELX:
02351 Vnm_print(2, "Vfetk_fillArray: can't write out x-shifted diel!\n");
02352 return 0;
02353 break;
02354
02355 case VDT_DIELY:
02356 Vnm_print(2, "Vfetk_fillArray: can't write out y-shifted diel!\n");
02357 return 0;
02358 break;
02359
02360 case VDT_DIELZ:
02361 Vnm_print(2, "Vfetk_fillArray: can't write out z-shifted diel!\n");
02362 return 0;
02363 break;
02364
02365 case VDT_KAPPA:
02366 Vnm_print(2, "Vfetk_fillArray: can't write out kappa!\n");
02367 return 0;
02368 break;
02369
02370 default:
02371 Vnm_print(2, "Vfetk_fillArray: invalid data type (%d)!\n", type);
02372 return 0;
02373 break;
02374 }
02375
02376 return 1;
02377 }
02378
02379 VPUBLIC int Vfetk_write(Vfetk *thee, const char *iodev, const char *iofmt,
02380 const char *thost, const char *fname, Bvec *vec, Vdata_Format format) {
02381
02382 int i, j, ichop;
02383 Aprx *aprx;
02384 Gem *gm;
02385 Vio *sock;
02386
02387 VASSERT(thee != VNULL);
02388 aprx = thee->aprx;
02389 gm = thee->gm;
02390
02391 sock = Vio_ctor(iodev,iofmt,thost,fname,"w");
02392 if (sock == VNULL) {
02393 Vnm_print(2, "Vfetk_write: Problem opening virtual socket %s\n",
02394 fname);
02395 return 0;
02396 }
02397 if (Vio_connect(sock, 0) < 0) {
02398 Vnm_print(2, "Vfetk_write: Problem connecting to virtual socket %s\n",
02399 fname);
02400 return 0;
02401 }
02402
02403
02404 if (Bvec_numRB(vec, 0) != Gem_numVV(gm)) {
02405 Vnm_print(2, "Vfetk_fillArray: insufficient space in Bvec!\n");
02406 Vnm_print(2, "Vfetk_fillArray: Have %d, need %d!\n", Bvec_numRB(vec, 0),
02407 Gem_numVV(gm));
02408 return 0;
02409 }
02410
02411 switch (format) {
02412
02413 case VDF_DX:
02414 Aprx_writeSOL(aprx, sock, vec, "DX");
02415 break;
02416 case VDF_AVS:
02417 Aprx_writeSOL(aprx, sock, vec, "UCD");
02418 break;
02419 case VDF_UHBD:
02420 Vnm_print(2, "Vfetk_write: UHBD format not supported!\n");
02421 return 0;
02422 default:
02423 Vnm_print(2, "Vfetk_write: Invalid data format (%d)!\n", format);
02424 return 0;
02425 }
02426
02427
02428 Vio_connectFree(sock);
02429 Vio_dtor(&sock);
02430
02431 return 1;
02432 }
02433
02434 #endif